|
Belos
Version of the Day
|
00001 //@HEADER 00002 // ************************************************************************ 00003 // 00004 // Belos: Block Linear Solvers Package 00005 // Copyright 2004 Sandia Corporation 00006 // 00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00008 // the U.S. Government retains certain rights in this software. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // ************************************************************************ 00040 //@HEADER 00041 00042 00047 #ifndef BELOS_DGKS_ORTHOMANAGER_HPP 00048 #define BELOS_DGKS_ORTHOMANAGER_HPP 00049 00057 // #define ORTHO_DEBUG 00058 00059 #include "BelosConfigDefs.hpp" 00060 #include "BelosMultiVecTraits.hpp" 00061 #include "BelosOperatorTraits.hpp" 00062 #include "BelosMatOrthoManager.hpp" 00063 00064 #include "Teuchos_as.hpp" 00065 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp" 00066 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00067 #include "Teuchos_TimeMonitor.hpp" 00068 #endif // BELOS_TEUCHOS_TIME_MONITOR 00069 00070 namespace Belos { 00071 00076 template<class ScalarType> 00077 Teuchos::RCP<const Teuchos::ParameterList> TEUCHOS_DEPRECATED 00078 getDefaultDgksParameters() 00079 { 00080 using Teuchos::as; 00081 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType magnitude_type; 00082 typedef Teuchos::ScalarTraits<magnitude_type> STM; 00083 00084 // Default parameter values for DGKS orthogonalization. 00085 // Documentation will be embedded in the parameter list. 00086 const int defaultMaxNumOrthogPasses = 2; 00087 const magnitude_type one = STM::one(); 00088 const magnitude_type eps = STM::eps(); 00089 const magnitude_type squareRootTwo = 00090 STM::squareroot (as<magnitude_type> (2)); 00091 const magnitude_type defaultBlkTol = 00092 as<magnitude_type> (10) * STM::squareroot (eps); 00093 const magnitude_type defaultDepTol = one / squareRootTwo; 00094 const magnitude_type defaultSingTol = as<magnitude_type> (10) * eps; 00095 00096 Teuchos::RCP<Teuchos::ParameterList> params = 00097 Teuchos::parameterList ("DGKS"); 00098 params->set ("maxNumOrthogPasses", defaultMaxNumOrthogPasses, 00099 "Maximum number of orthogonalization passes " 00100 "(includes the first). Default is 2, since " 00101 "\"twice is enough\" for Krylov methods."); 00102 params->set ("blkTol", defaultBlkTol, 00103 "Block reorthogonalization threshhold."); 00104 params->set ("depTol", defaultDepTol, 00105 "(Non-block) reorthogonalization threshold."); 00106 params->set ("singTol", defaultSingTol, 00107 "Singular block detection threshold."); 00108 return params; 00109 } 00110 00116 template<class ScalarType> 00117 Teuchos::RCP<const Teuchos::ParameterList> TEUCHOS_DEPRECATED 00118 getFastDgksParameters() 00119 { 00120 using Teuchos::ParameterList; 00121 using Teuchos::RCP; 00122 using Teuchos::rcp; 00123 using Teuchos::ScalarTraits; 00124 typedef typename ScalarTraits<ScalarType>::magnitudeType magnitude_type; 00125 typedef ScalarTraits<magnitude_type> STM; 00126 00127 RCP<const ParameterList> defaultParams = getDefaultDgksParameters<ScalarType> (); 00128 // Start with a clone of the default parameters. 00129 RCP<ParameterList> params = rcp (new ParameterList (*defaultParams)); 00130 00131 const int maxBlkOrtho = 1; 00132 const magnitude_type blkTol = STM::zero(); 00133 const magnitude_type depTol = STM::zero(); 00134 const magnitude_type singTol = STM::zero(); 00135 00136 params->set ("maxNumOrthogPasses", maxBlkOrtho); 00137 params->set ("blkTol", blkTol); 00138 params->set ("depTol", depTol); 00139 params->set ("singTol", singTol); 00140 00141 return params; 00142 } 00143 00153 template<class ScalarType> 00154 void TEUCHOS_DEPRECATED 00155 readDgksParameters (const Teuchos::RCP<const Teuchos::ParameterList>& params, 00156 int& maxNumOrthogPasses, 00157 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType& blkTol, 00158 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType& depTol, 00159 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType& singTol) 00160 { 00161 using Teuchos::ParameterList; 00162 using Teuchos::RCP; 00163 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType magnitude_type; 00164 const magnitude_type zero = Teuchos::ScalarTraits<magnitude_type>::zero(); 00165 RCP<const ParameterList> defaultParams = getDefaultDgksParameters<ScalarType>(); 00166 00167 // Using temporary variables and fetching all values before 00168 // setting the output arguments ensures the strong exception 00169 // guarantee for this function: if an exception is thrown, no 00170 // externally visible side effects (in this case, setting the 00171 // output arguments) have taken place. 00172 int _maxNumOrthogPasses; 00173 magnitude_type _blkTol, _depTol, _singTol; 00174 if (params.is_null()) 00175 { 00176 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00177 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00178 _depTol = defaultParams->get<magnitude_type> ("depTol"); 00179 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00180 } 00181 else 00182 { 00183 try { 00184 _maxNumOrthogPasses = params->get<int> ("maxNumOrthogPasses"); 00185 if (_maxNumOrthogPasses < 1) 00186 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00187 } catch (Teuchos::Exceptions::InvalidParameter&) { 00188 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00189 } 00190 00191 try { 00192 _blkTol = params->get<magnitude_type> ("blkTol"); 00193 if (_blkTol < zero) 00194 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00195 } catch (Teuchos::Exceptions::InvalidParameter&) { 00196 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00197 } 00198 00199 try { 00200 _depTol = params->get<magnitude_type> ("depTol"); 00201 if (_depTol < zero) 00202 _depTol = defaultParams->get<magnitude_type> ("depTol"); 00203 } catch (Teuchos::Exceptions::InvalidParameter&) { 00204 _depTol = defaultParams->get<magnitude_type> ("depTol"); 00205 } 00206 00207 try { 00208 _singTol = params->get<magnitude_type> ("singTol"); 00209 if (_singTol < zero) 00210 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00211 } catch (Teuchos::Exceptions::InvalidParameter&) { 00212 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00213 } 00214 } 00215 maxNumOrthogPasses = _maxNumOrthogPasses; 00216 blkTol = _blkTol; 00217 depTol = _depTol; 00218 singTol = _singTol; 00219 } 00220 00221 template<class ScalarType, class MV, class OP> 00222 class DGKSOrthoManager : 00223 public MatOrthoManager<ScalarType,MV,OP>, 00224 public Teuchos::ParameterListAcceptorDefaultBase 00225 { 00226 private: 00227 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType; 00228 typedef typename Teuchos::ScalarTraits<MagnitudeType> MGT; 00229 typedef Teuchos::ScalarTraits<ScalarType> SCT; 00230 typedef MultiVecTraits<ScalarType,MV> MVT; 00231 typedef MultiVecTraitsExt<ScalarType,MV> MVText; 00232 typedef OperatorTraits<ScalarType,MV,OP> OPT; 00233 00234 public: 00236 00237 00239 DGKSOrthoManager( const std::string& label = "Belos", 00240 Teuchos::RCP<const OP> Op = Teuchos::null, 00241 const int max_blk_ortho = 2, 00242 const MagnitudeType blk_tol = 10*MGT::squareroot( MGT::eps() ), 00243 const MagnitudeType dep_tol = MGT::one()/MGT::squareroot( 2*MGT::one() ), 00244 const MagnitudeType sing_tol = 10*MGT::eps() ) 00245 : MatOrthoManager<ScalarType,MV,OP>(Op), 00246 max_blk_ortho_( max_blk_ortho ), 00247 blk_tol_( blk_tol ), 00248 dep_tol_( dep_tol ), 00249 sing_tol_( sing_tol ), 00250 label_( label ) 00251 { 00252 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00253 std::string orthoLabel = label_ + ": Orthogonalization"; 00254 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter( orthoLabel ); 00255 #endif 00256 } 00257 00259 DGKSOrthoManager (const Teuchos::RCP<Teuchos::ParameterList>& plist, 00260 const std::string& label = "Belos", 00261 Teuchos::RCP<const OP> Op = Teuchos::null) 00262 : MatOrthoManager<ScalarType,MV,OP>(Op), 00263 max_blk_ortho_ (2), 00264 blk_tol_ (Teuchos::as<MagnitudeType>(10) * MGT::squareroot(MGT::eps())), 00265 dep_tol_ (MGT::one() / MGT::squareroot (Teuchos::as<MagnitudeType>(2))), 00266 sing_tol_ (Teuchos::as<MagnitudeType>(10) * MGT::eps()), 00267 label_( label ) 00268 { 00269 setParameterList (plist); 00270 00271 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00272 std::string orthoLabel = label_ + ": Orthogonalization"; 00273 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter( orthoLabel ); 00274 #endif 00275 } 00276 00278 ~DGKSOrthoManager() {} 00280 00282 00283 00284 void 00285 setParameterList (const Teuchos::RCP<Teuchos::ParameterList>& plist) 00286 { 00287 using Teuchos::ParameterList; 00288 using Teuchos::parameterList; 00289 using Teuchos::RCP; 00290 00291 RCP<const ParameterList> defaultParams = getValidParameters(); 00292 RCP<ParameterList> params; 00293 if (plist.is_null()) { 00294 // No need to validate default parameters. 00295 params = parameterList (*defaultParams); 00296 } else { 00297 params = plist; 00298 params->validateParametersAndSetDefaults (*defaultParams); 00299 } 00300 00301 // Using temporary variables and fetching all values before 00302 // setting the output arguments ensures the strong exception 00303 // guarantee for this function: if an exception is thrown, no 00304 // externally visible side effects (in this case, setting the 00305 // output arguments) have taken place. 00306 const int maxNumOrthogPasses = params->get<int> ("maxNumOrthogPasses"); 00307 const MagnitudeType blkTol = params->get<MagnitudeType> ("blkTol"); 00308 const MagnitudeType depTol = params->get<MagnitudeType> ("depTol"); 00309 const MagnitudeType singTol = params->get<MagnitudeType> ("singTol"); 00310 00311 max_blk_ortho_ = maxNumOrthogPasses; 00312 blk_tol_ = blkTol; 00313 dep_tol_ = depTol; 00314 sing_tol_ = singTol; 00315 00316 setMyParamList (params); 00317 } 00318 00319 Teuchos::RCP<const Teuchos::ParameterList> 00320 getValidParameters () const 00321 { 00322 using Teuchos::as; 00323 using Teuchos::ParameterList; 00324 using Teuchos::parameterList; 00325 using Teuchos::RCP; 00326 00327 if (defaultParams_.is_null()) { 00328 RCP<ParameterList> params = parameterList ("DGKS"); 00329 const MagnitudeType eps = MGT::eps (); 00330 00331 // Default parameter values for DGKS orthogonalization. 00332 // Documentation will be embedded in the parameter list. 00333 const int defaultMaxNumOrthogPasses = 2; 00334 const MagnitudeType defaultBlkTol = 00335 as<MagnitudeType> (10) * MGT::squareroot (eps); 00336 const MagnitudeType defaultDepTol = 00337 MGT::one() / MGT::squareroot (as<MagnitudeType> (2)); 00338 00339 const MagnitudeType defaultSingTol = as<MagnitudeType> (10) * eps; 00340 00341 params->set ("maxNumOrthogPasses", defaultMaxNumOrthogPasses, 00342 "Maximum number of orthogonalization passes (includes the " 00343 "first). Default is 2, since \"twice is enough\" for Krylov " 00344 "methods."); 00345 params->set ("blkTol", defaultBlkTol, "Block reorthogonalization " 00346 "threshhold."); 00347 params->set ("depTol", defaultDepTol, 00348 "(Non-block) reorthogonalization threshold."); 00349 params->set ("singTol", defaultSingTol, "Singular block detection " 00350 "threshold."); 00351 defaultParams_ = params; 00352 } 00353 return defaultParams_; 00354 } 00355 00357 00358 Teuchos::RCP<const Teuchos::ParameterList> 00359 getFastParameters() const 00360 { 00361 using Teuchos::ParameterList; 00362 using Teuchos::RCP; 00363 using Teuchos::rcp; 00364 00365 RCP<const ParameterList> defaultParams = getValidParameters (); 00366 // Start with a clone of the default parameters. 00367 RCP<ParameterList> params = rcp (new ParameterList (*defaultParams)); 00368 00369 const int maxBlkOrtho = 1; 00370 const MagnitudeType blkTol = MGT::zero(); 00371 const MagnitudeType depTol = MGT::zero(); 00372 const MagnitudeType singTol = MGT::zero(); 00373 00374 params->set ("maxNumOrthogPasses", maxBlkOrtho); 00375 params->set ("blkTol", blkTol); 00376 params->set ("depTol", depTol); 00377 params->set ("singTol", singTol); 00378 00379 return params; 00380 } 00381 00383 00384 00386 void setBlkTol( const MagnitudeType blk_tol ) { 00387 // Update the parameter list as well. 00388 Teuchos::RCP<Teuchos::ParameterList> params = getNonconstParameterList(); 00389 if (! params.is_null()) { 00390 // If it's null, then we haven't called setParameterList() 00391 // yet. It's entirely possible to construct the parameter 00392 // list on demand, so we don't try to create the parameter 00393 // list here. 00394 params->set ("blkTol", blk_tol); 00395 } 00396 blk_tol_ = blk_tol; 00397 } 00398 00400 void setDepTol( const MagnitudeType dep_tol ) { 00401 // Update the parameter list as well. 00402 Teuchos::RCP<Teuchos::ParameterList> params = getNonconstParameterList(); 00403 if (! params.is_null()) { 00404 params->set ("depTol", dep_tol); 00405 } 00406 dep_tol_ = dep_tol; 00407 } 00408 00410 void setSingTol( const MagnitudeType sing_tol ) { 00411 // Update the parameter list as well. 00412 Teuchos::RCP<Teuchos::ParameterList> params = getNonconstParameterList(); 00413 if (! params.is_null()) { 00414 params->set ("singTol", sing_tol); 00415 } 00416 sing_tol_ = sing_tol; 00417 } 00418 00420 MagnitudeType getBlkTol() const { return blk_tol_; } 00421 00423 MagnitudeType getDepTol() const { return dep_tol_; } 00424 00426 MagnitudeType getSingTol() const { return sing_tol_; } 00427 00429 00430 00432 00433 00461 void project ( MV &X, Teuchos::RCP<MV> MX, 00462 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00463 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00464 00465 00468 void project ( MV &X, 00469 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00470 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const { 00471 project(X,Teuchos::null,C,Q); 00472 } 00473 00474 00475 00500 int normalize ( MV &X, Teuchos::RCP<MV> MX, 00501 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B) const; 00502 00503 00506 int normalize ( MV &X, Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const { 00507 return normalize(X,Teuchos::null,B); 00508 } 00509 00510 protected: 00511 00567 virtual int 00568 projectAndNormalizeWithMxImpl (MV &X, 00569 Teuchos::RCP<MV> MX, 00570 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00571 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00572 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00573 00574 public: 00576 00577 00578 00584 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00585 orthonormError(const MV &X) const { 00586 return orthonormError(X,Teuchos::null); 00587 } 00588 00595 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00596 orthonormError(const MV &X, Teuchos::RCP<const MV> MX) const; 00597 00601 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00602 orthogError(const MV &X1, const MV &X2) const { 00603 return orthogError(X1,Teuchos::null,X2); 00604 } 00605 00610 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00611 orthogError(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const; 00612 00614 00616 00617 00620 void setLabel(const std::string& label); 00621 00624 const std::string& getLabel() const { return label_; } 00625 00627 00628 private: 00629 00631 int max_blk_ortho_; 00633 MagnitudeType blk_tol_; 00635 MagnitudeType dep_tol_; 00637 MagnitudeType sing_tol_; 00638 00640 std::string label_; 00641 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00642 Teuchos::RCP<Teuchos::Time> timerOrtho_; 00643 #endif // BELOS_TEUCHOS_TIME_MONITOR 00644 00646 mutable Teuchos::RCP<Teuchos::ParameterList> defaultParams_; 00647 00649 int findBasis(MV &X, Teuchos::RCP<MV> MX, 00650 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > C, 00651 bool completeBasis, int howMany = -1 ) const; 00652 00654 bool blkOrtho ( MV &X, Teuchos::RCP<MV> MX, 00655 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00656 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00657 00671 int blkOrthoSing ( MV &X, Teuchos::RCP<MV> MX, 00672 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00673 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00674 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ) const; 00675 }; 00676 00678 // Set the label for this orthogonalization manager and create new timers if it's changed 00679 template<class ScalarType, class MV, class OP> 00680 void DGKSOrthoManager<ScalarType,MV,OP>::setLabel(const std::string& label) 00681 { 00682 if (label != label_) { 00683 label_ = label; 00684 std::string orthoLabel = label_ + ": Orthogonalization"; 00685 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00686 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel); 00687 #endif 00688 } 00689 } 00690 00692 // Compute the distance from orthonormality 00693 template<class ScalarType, class MV, class OP> 00694 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00695 DGKSOrthoManager<ScalarType,MV,OP>::orthonormError(const MV &X, Teuchos::RCP<const MV> MX) const { 00696 const ScalarType ONE = SCT::one(); 00697 int rank = MVT::GetNumberVecs(X); 00698 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank); 00699 MatOrthoManager<ScalarType,MV,OP>::innerProd(X,X,MX,xTx); 00700 for (int i=0; i<rank; i++) { 00701 xTx(i,i) -= ONE; 00702 } 00703 return xTx.normFrobenius(); 00704 } 00705 00707 // Compute the distance from orthogonality 00708 template<class ScalarType, class MV, class OP> 00709 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00710 DGKSOrthoManager<ScalarType,MV,OP>::orthogError(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const { 00711 int r1 = MVT::GetNumberVecs(X1); 00712 int r2 = MVT::GetNumberVecs(X2); 00713 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r2,r1); 00714 MatOrthoManager<ScalarType,MV,OP>::innerProd(X2,X1,MX1,xTx); 00715 return xTx.normFrobenius(); 00716 } 00717 00719 // Find an Op-orthonormal basis for span(X) - span(W) 00720 template<class ScalarType, class MV, class OP> 00721 int 00722 DGKSOrthoManager<ScalarType, MV, OP>:: 00723 projectAndNormalizeWithMxImpl (MV &X, 00724 Teuchos::RCP<MV> MX, 00725 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00726 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00727 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const 00728 { 00729 using Teuchos::Array; 00730 using Teuchos::null; 00731 using Teuchos::is_null; 00732 using Teuchos::RCP; 00733 using Teuchos::rcp; 00734 using Teuchos::SerialDenseMatrix; 00735 typedef SerialDenseMatrix< int, ScalarType > serial_dense_matrix_type; 00736 typedef typename Array< RCP< const MV > >::size_type size_type; 00737 00738 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00739 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00740 #endif 00741 00742 ScalarType ONE = SCT::one(); 00743 ScalarType ZERO = SCT::zero(); 00744 00745 int nq = Q.size(); 00746 int xc = MVT::GetNumberVecs( X ); 00747 ptrdiff_t xr = MVText::GetGlobalLength( X ); 00748 int rank = xc; 00749 00750 // If the user doesn't want to store the normalization 00751 // coefficients, allocate some local memory for them. This will 00752 // go away at the end of this method. 00753 if (is_null (B)) { 00754 B = rcp (new serial_dense_matrix_type (xc, xc)); 00755 } 00756 // Likewise, if the user doesn't want to store the projection 00757 // coefficients, allocate some local memory for them. Also make 00758 // sure that all the entries of C are the right size. We're going 00759 // to overwrite them anyway, so we don't have to worry about the 00760 // contents (other than to resize them if they are the wrong 00761 // size). 00762 if (C.size() < nq) 00763 C.resize (nq); 00764 for (size_type k = 0; k < nq; ++k) 00765 { 00766 const int numRows = MVT::GetNumberVecs (*Q[k]); 00767 const int numCols = xc; // Number of vectors in X 00768 00769 if (is_null (C[k])) 00770 C[k] = rcp (new serial_dense_matrix_type (numRows, numCols)); 00771 else if (C[k]->numRows() != numRows || C[k]->numCols() != numCols) 00772 { 00773 int err = C[k]->reshape (numRows, numCols); 00774 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, 00775 "DGKS orthogonalization: failed to reshape " 00776 "C[" << k << "] (the array of block " 00777 "coefficients resulting from projecting X " 00778 "against Q[1:" << nq << "])."); 00779 } 00780 } 00781 00782 /****** DO NO MODIFY *MX IF _hasOp == false ******/ 00783 if (this->_hasOp) { 00784 if (MX == Teuchos::null) { 00785 // we need to allocate space for MX 00786 MX = MVT::Clone(X,MVT::GetNumberVecs(X)); 00787 OPT::Apply(*(this->_Op),X,*MX); 00788 } 00789 } 00790 else { 00791 // Op == I --> MX = X (ignore it if the user passed it in) 00792 MX = Teuchos::rcp( &X, false ); 00793 } 00794 00795 int mxc = MVT::GetNumberVecs( *MX ); 00796 ptrdiff_t mxr = MVText::GetGlobalLength( *MX ); 00797 00798 // short-circuit 00799 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument, "Belos::DGKSOrthoManager::projectAndNormalize(): X must be non-empty" ); 00800 00801 int numbas = 0; 00802 for (int i=0; i<nq; i++) { 00803 numbas += MVT::GetNumberVecs( *Q[i] ); 00804 } 00805 00806 // check size of B 00807 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument, 00808 "Belos::DGKSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of B" ); 00809 // check size of X and MX 00810 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument, 00811 "Belos::DGKSOrthoManager::projectAndNormalize(): MVT returned negative dimensions for X,MX" ); 00812 // check size of X w.r.t. MX 00813 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument, 00814 "Belos::DGKSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of MX" ); 00815 // check feasibility 00816 //TEUCHOS_TEST_FOR_EXCEPTION( numbas+xc > xr, std::invalid_argument, 00817 // "Belos::DGKSOrthoManager::projectAndNormalize(): Orthogonality constraints not feasible" ); 00818 00819 // Some flags for checking dependency returns from the internal orthogonalization methods 00820 bool dep_flg = false; 00821 00822 // Make a temporary copy of X and MX, just in case a block dependency is detected. 00823 Teuchos::RCP<MV> tmpX, tmpMX; 00824 tmpX = MVT::CloneCopy(X); 00825 if (this->_hasOp) { 00826 tmpMX = MVT::CloneCopy(*MX); 00827 } 00828 00829 // Use the cheaper block orthogonalization. 00830 dep_flg = blkOrtho( X, MX, C, Q ); 00831 00832 // If a dependency has been detected in this block, then perform 00833 // the more expensive single-vector orthogonalization. 00834 if (dep_flg) { 00835 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q ); 00836 00837 // Copy tmpX back into X. 00838 MVT::MvAddMv( ONE, *tmpX, ZERO, *tmpX, X ); 00839 if (this->_hasOp) { 00840 MVT::MvAddMv( ONE, *tmpMX, ZERO, *tmpMX, *MX ); 00841 } 00842 } 00843 else { 00844 // There is no dependency, so orthonormalize new block X 00845 rank = findBasis( X, MX, B, false ); 00846 if (rank < xc) { 00847 // A dependency was found during orthonormalization of X, 00848 // rerun orthogonalization using more expensive single-vector orthogonalization. 00849 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q ); 00850 00851 // Copy tmpX back into X. 00852 MVT::MvAddMv( ONE, *tmpX, ZERO, *tmpX, X ); 00853 if (this->_hasOp) { 00854 MVT::MvAddMv( ONE, *tmpMX, ZERO, *tmpMX, *MX ); 00855 } 00856 } 00857 } 00858 00859 // this should not raise an std::exception; but our post-conditions oblige us to check 00860 TEUCHOS_TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error, 00861 "Belos::DGKSOrthoManager::projectAndNormalize(): Debug error in rank variable." ); 00862 00863 // Return the rank of X. 00864 return rank; 00865 } 00866 00867 00868 00870 // Find an Op-orthonormal basis for span(X), with rank numvectors(X) 00871 template<class ScalarType, class MV, class OP> 00872 int DGKSOrthoManager<ScalarType, MV, OP>::normalize( 00873 MV &X, Teuchos::RCP<MV> MX, 00874 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const { 00875 00876 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00877 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00878 #endif 00879 00880 // call findBasis, with the instruction to try to generate a basis of rank numvecs(X) 00881 return findBasis(X, MX, B, true); 00882 } 00883 00884 00885 00887 template<class ScalarType, class MV, class OP> 00888 void DGKSOrthoManager<ScalarType, MV, OP>::project( 00889 MV &X, Teuchos::RCP<MV> MX, 00890 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00891 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const { 00892 // For the inner product defined by the operator Op or the identity (Op == 0) 00893 // -> Orthogonalize X against each Q[i] 00894 // Modify MX accordingly 00895 // 00896 // Note that when Op is 0, MX is not referenced 00897 // 00898 // Parameter variables 00899 // 00900 // X : Vectors to be transformed 00901 // 00902 // MX : Image of the block vector X by the mass matrix 00903 // 00904 // Q : Bases to orthogonalize against. These are assumed orthonormal, mutually and independently. 00905 // 00906 00907 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00908 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00909 #endif 00910 00911 int xc = MVT::GetNumberVecs( X ); 00912 ptrdiff_t xr = MVText::GetGlobalLength( X ); 00913 int nq = Q.size(); 00914 std::vector<int> qcs(nq); 00915 // short-circuit 00916 if (nq == 0 || xc == 0 || xr == 0) { 00917 return; 00918 } 00919 ptrdiff_t qr = MVText::GetGlobalLength ( *Q[0] ); 00920 // if we don't have enough C, expand it with null references 00921 // if we have too many, resize to throw away the latter ones 00922 // if we have exactly as many as we have Q, this call has no effect 00923 C.resize(nq); 00924 00925 00926 /****** DO NO MODIFY *MX IF _hasOp == false ******/ 00927 if (this->_hasOp) { 00928 if (MX == Teuchos::null) { 00929 // we need to allocate space for MX 00930 MX = MVT::Clone(X,MVT::GetNumberVecs(X)); 00931 OPT::Apply(*(this->_Op),X,*MX); 00932 } 00933 } 00934 else { 00935 // Op == I --> MX = X (ignore it if the user passed it in) 00936 MX = Teuchos::rcp( &X, false ); 00937 } 00938 int mxc = MVT::GetNumberVecs( *MX ); 00939 ptrdiff_t mxr = MVText::GetGlobalLength( *MX ); 00940 00941 // check size of X and Q w.r.t. common sense 00942 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument, 00943 "Belos::DGKSOrthoManager::project(): MVT returned negative dimensions for X,MX" ); 00944 // check size of X w.r.t. MX and Q 00945 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr || xr!=qr, std::invalid_argument, 00946 "Belos::DGKSOrthoManager::project(): Size of X not consistant with MX,Q" ); 00947 00948 // tally up size of all Q and check/allocate C 00949 int baslen = 0; 00950 for (int i=0; i<nq; i++) { 00951 TEUCHOS_TEST_FOR_EXCEPTION( MVText::GetGlobalLength( *Q[i] ) != qr, std::invalid_argument, 00952 "Belos::DGKSOrthoManager::project(): Q lengths not mutually consistant" ); 00953 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 00954 TEUCHOS_TEST_FOR_EXCEPTION( qr < qcs[i], std::invalid_argument, 00955 "Belos::DGKSOrthoManager::project(): Q has less rows than columns" ); 00956 baslen += qcs[i]; 00957 00958 // check size of C[i] 00959 if ( C[i] == Teuchos::null ) { 00960 C[i] = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) ); 00961 } 00962 else { 00963 TEUCHOS_TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc , std::invalid_argument, 00964 "Belos::DGKSOrthoManager::project(): Size of Q not consistant with size of C" ); 00965 } 00966 } 00967 00968 // Use the cheaper block orthogonalization, don't check for rank deficiency. 00969 blkOrtho( X, MX, C, Q ); 00970 00971 } 00972 00974 // Find an Op-orthonormal basis for span(X), with the option of extending the subspace so that 00975 // the rank is numvectors(X) 00976 template<class ScalarType, class MV, class OP> 00977 int DGKSOrthoManager<ScalarType, MV, OP>::findBasis( 00978 MV &X, Teuchos::RCP<MV> MX, 00979 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00980 bool completeBasis, int howMany ) const { 00981 // For the inner product defined by the operator Op or the identity (Op == 0) 00982 // -> Orthonormalize X 00983 // Modify MX accordingly 00984 // 00985 // Note that when Op is 0, MX is not referenced 00986 // 00987 // Parameter variables 00988 // 00989 // X : Vectors to be orthonormalized 00990 // 00991 // MX : Image of the multivector X under the operator Op 00992 // 00993 // Op : Pointer to the operator for the inner product 00994 // 00995 // 00996 00997 const ScalarType ONE = SCT::one(); 00998 const MagnitudeType ZERO = SCT::magnitude(SCT::zero()); 00999 01000 int xc = MVT::GetNumberVecs( X ); 01001 ptrdiff_t xr = MVText::GetGlobalLength( X ); 01002 01003 if (howMany == -1) { 01004 howMany = xc; 01005 } 01006 01007 /******************************************************* 01008 * If _hasOp == false, we will not reference MX below * 01009 *******************************************************/ 01010 01011 // if Op==null, MX == X (via pointer) 01012 // Otherwise, either the user passed in MX or we will allocated and compute it 01013 if (this->_hasOp) { 01014 if (MX == Teuchos::null) { 01015 // we need to allocate space for MX 01016 MX = MVT::Clone(X,xc); 01017 OPT::Apply(*(this->_Op),X,*MX); 01018 } 01019 } 01020 01021 /* if the user doesn't want to store the coefficienets, 01022 * allocate some local memory for them 01023 */ 01024 if ( B == Teuchos::null ) { 01025 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) ); 01026 } 01027 01028 int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc; 01029 ptrdiff_t mxr = (this->_hasOp) ? MVText::GetGlobalLength( *MX ) : xr; 01030 01031 // check size of C, B 01032 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument, 01033 "Belos::DGKSOrthoManager::findBasis(): X must be non-empty" ); 01034 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument, 01035 "Belos::DGKSOrthoManager::findBasis(): Size of X not consistant with size of B" ); 01036 TEUCHOS_TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument, 01037 "Belos::DGKSOrthoManager::findBasis(): Size of X not consistant with size of MX" ); 01038 TEUCHOS_TEST_FOR_EXCEPTION( static_cast<ptrdiff_t>(xc) > xr, std::invalid_argument, 01039 "Belos::DGKSOrthoManager::findBasis(): Size of X not feasible for normalization" ); 01040 TEUCHOS_TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::invalid_argument, 01041 "Belos::DGKSOrthoManager::findBasis(): Invalid howMany parameter" ); 01042 01043 /* xstart is which column we are starting the process with, based on howMany 01044 * columns before xstart are assumed to be Op-orthonormal already 01045 */ 01046 int xstart = xc - howMany; 01047 01048 for (int j = xstart; j < xc; j++) { 01049 01050 // numX is 01051 // * number of currently orthonormal columns of X 01052 // * the index of the current column of X 01053 int numX = j; 01054 bool addVec = false; 01055 01056 // Get a view of the vector currently being worked on. 01057 std::vector<int> index(1); 01058 index[0] = numX; 01059 Teuchos::RCP<MV> Xj = MVT::CloneViewNonConst( X, index ); 01060 Teuchos::RCP<MV> MXj; 01061 if ((this->_hasOp)) { 01062 // MXj is a view of the current vector in MX 01063 MXj = MVT::CloneViewNonConst( *MX, index ); 01064 } 01065 else { 01066 // MXj is a pointer to Xj, and MUST NOT be modified 01067 MXj = Xj; 01068 } 01069 01070 // Get a view of the previous vectors. 01071 std::vector<int> prev_idx( numX ); 01072 Teuchos::RCP<const MV> prevX, prevMX; 01073 01074 if (numX > 0) { 01075 for (int i=0; i<numX; i++) { 01076 prev_idx[i] = i; 01077 } 01078 prevX = MVT::CloneView( X, prev_idx ); 01079 if (this->_hasOp) { 01080 prevMX = MVT::CloneView( *MX, prev_idx ); 01081 } 01082 } 01083 01084 // Make storage for these Gram-Schmidt iterations. 01085 Teuchos::SerialDenseMatrix<int,ScalarType> product(numX, 1); 01086 std::vector<ScalarType> oldDot( 1 ), newDot( 1 ); 01087 // 01088 // Save old MXj vector and compute Op-norm 01089 // 01090 Teuchos::RCP<MV> oldMXj = MVT::CloneCopy( *MXj ); 01091 MVT::MvDot( *Xj, *MXj, oldDot ); 01092 // Xj^H Op Xj should be real and positive, by the hermitian positive definiteness of Op 01093 TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(oldDot[0]) < ZERO, OrthoError, 01094 "Belos::DGKSOrthoManager::findBasis(): Negative definiteness discovered in inner product" ); 01095 01096 if (numX > 0) { 01097 // Apply the first step of Gram-Schmidt 01098 01099 // product <- prevX^T MXj 01100 MatOrthoManager<ScalarType,MV,OP>::innerProd(*prevX,*Xj,MXj,product); 01101 01102 // Xj <- Xj - prevX prevX^T MXj 01103 // = Xj - prevX product 01104 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *Xj ); 01105 01106 // Update MXj 01107 if (this->_hasOp) { 01108 // MXj <- Op*Xj_new 01109 // = Op*(Xj_old - prevX prevX^T MXj) 01110 // = MXj - prevMX product 01111 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *MXj ); 01112 } 01113 01114 // Compute new Op-norm 01115 MVT::MvDot( *Xj, *MXj, newDot ); 01116 01117 // Check if a correction is needed. 01118 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(dep_tol_*oldDot[0]) ) { 01119 // Apply the second step of Gram-Schmidt 01120 // This is the same as above 01121 Teuchos::SerialDenseMatrix<int,ScalarType> P2(numX,1); 01122 01123 MatOrthoManager<ScalarType,MV,OP>::innerProd(*prevX,*Xj,MXj,P2); 01124 product += P2; 01125 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj ); 01126 if ((this->_hasOp)) { 01127 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj ); 01128 } 01129 } // if (newDot[0] < dep_tol_*oldDot[0]) 01130 01131 } // if (numX > 0) 01132 01133 // Compute Op-norm with old MXj 01134 MVT::MvDot( *Xj, *oldMXj, newDot ); 01135 01136 // Check to see if the new vector is dependent. 01137 if (completeBasis) { 01138 // 01139 // We need a complete basis, so add random vectors if necessary 01140 // 01141 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(sing_tol_*oldDot[0]) ) { 01142 01143 // Add a random vector and orthogonalize it against previous vectors in block. 01144 addVec = true; 01145 #ifdef ORTHO_DEBUG 01146 std::cout << "Belos::DGKSOrthoManager::findBasis() --> Random for column " << numX << std::endl; 01147 #endif 01148 // 01149 Teuchos::RCP<MV> tempXj = MVT::Clone( X, 1 ); 01150 Teuchos::RCP<MV> tempMXj; 01151 MVT::MvRandom( *tempXj ); 01152 if (this->_hasOp) { 01153 tempMXj = MVT::Clone( X, 1 ); 01154 OPT::Apply( *(this->_Op), *tempXj, *tempMXj ); 01155 } 01156 else { 01157 tempMXj = tempXj; 01158 } 01159 MVT::MvDot( *tempXj, *tempMXj, oldDot ); 01160 // 01161 for (int num_orth=0; num_orth<max_blk_ortho_; num_orth++){ 01162 MatOrthoManager<ScalarType,MV,OP>::innerProd(*prevX,*tempXj,tempMXj,product); 01163 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *tempXj ); 01164 if (this->_hasOp) { 01165 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *tempMXj ); 01166 } 01167 } 01168 // Compute new Op-norm 01169 MVT::MvDot( *tempXj, *tempMXj, newDot ); 01170 // 01171 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) { 01172 // Copy vector into current column of _basisvecs 01173 MVT::MvAddMv( ONE, *tempXj, ZERO, *tempXj, *Xj ); 01174 if (this->_hasOp) { 01175 MVT::MvAddMv( ONE, *tempMXj, ZERO, *tempMXj, *MXj ); 01176 } 01177 } 01178 else { 01179 return numX; 01180 } 01181 } 01182 } 01183 else { 01184 // 01185 // We only need to detect dependencies. 01186 // 01187 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*blk_tol_) ) { 01188 return numX; 01189 } 01190 } 01191 01192 // If we haven't left this method yet, then we can normalize the new vector Xj. 01193 // Normalize Xj. 01194 // Xj <- Xj / std::sqrt(newDot) 01195 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01196 01197 if (SCT::magnitude(diag) > ZERO) { 01198 MVT::MvAddMv( ONE/diag, *Xj, ZERO, *Xj, *Xj ); 01199 if (this->_hasOp) { 01200 // Update MXj. 01201 MVT::MvAddMv( ONE/diag, *MXj, ZERO, *MXj, *MXj ); 01202 } 01203 } 01204 01205 // If we've added a random vector, enter a zero in the j'th diagonal element. 01206 if (addVec) { 01207 (*B)(j,j) = ZERO; 01208 } 01209 else { 01210 (*B)(j,j) = diag; 01211 } 01212 01213 // Save the coefficients, if we are working on the original vector and not a randomly generated one 01214 if (!addVec) { 01215 for (int i=0; i<numX; i++) { 01216 (*B)(i,j) = product(i,0); 01217 } 01218 } 01219 01220 } // for (j = 0; j < xc; ++j) 01221 01222 return xc; 01223 } 01224 01226 // Routine to compute the block orthogonalization 01227 template<class ScalarType, class MV, class OP> 01228 bool 01229 DGKSOrthoManager<ScalarType, MV, OP>::blkOrtho ( MV &X, Teuchos::RCP<MV> MX, 01230 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 01231 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const 01232 { 01233 int nq = Q.size(); 01234 int xc = MVT::GetNumberVecs( X ); 01235 bool dep_flg = false; 01236 const ScalarType ONE = SCT::one(); 01237 01238 std::vector<int> qcs( nq ); 01239 for (int i=0; i<nq; i++) { 01240 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01241 } 01242 01243 // Perform the Gram-Schmidt transformation for a block of vectors 01244 01245 // Compute the initial Op-norms 01246 std::vector<ScalarType> oldDot( xc ); 01247 MVT::MvDot( X, *MX, oldDot ); 01248 01249 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq); 01250 // Define the product Q^T * (Op*X) 01251 for (int i=0; i<nq; i++) { 01252 // Multiply Q' with MX 01253 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,*C[i]); 01254 // Multiply by Q and subtract the result in X 01255 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X ); 01256 01257 // Update MX, with the least number of applications of Op as possible 01258 if (this->_hasOp) { 01259 if (xc <= qcs[i]) { 01260 OPT::Apply( *(this->_Op), X, *MX); 01261 } 01262 else { 01263 // this will possibly be used again below; don't delete it 01264 MQ[i] = MVT::Clone( *Q[i], qcs[i] ); 01265 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] ); 01266 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX ); 01267 } 01268 } 01269 } 01270 01271 // Do as many steps of classical Gram-Schmidt as required by max_blk_ortho_ 01272 for (int j = 1; j < max_blk_ortho_; ++j) { 01273 01274 for (int i=0; i<nq; i++) { 01275 Teuchos::SerialDenseMatrix<int,ScalarType> C2(*C[i]); 01276 01277 // Apply another step of classical Gram-Schmidt 01278 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,C2); 01279 *C[i] += C2; 01280 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X ); 01281 01282 // Update MX, with the least number of applications of Op as possible 01283 if (this->_hasOp) { 01284 if (MQ[i].get()) { 01285 // MQ was allocated and computed above; use it 01286 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX ); 01287 } 01288 else if (xc <= qcs[i]) { 01289 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01290 OPT::Apply( *(this->_Op), X, *MX); 01291 } 01292 } 01293 } // for (int i=0; i<nq; i++) 01294 } // for (int j = 0; j < max_blk_ortho; ++j) 01295 01296 // Compute new Op-norms 01297 std::vector<ScalarType> newDot(xc); 01298 MVT::MvDot( X, *MX, newDot ); 01299 01300 // Check to make sure the new block of vectors are not dependent on previous vectors 01301 for (int i=0; i<xc; i++){ 01302 if (SCT::magnitude(newDot[i]) < SCT::magnitude(oldDot[i] * blk_tol_)) { 01303 dep_flg = true; 01304 break; 01305 } 01306 } // end for (i=0;...) 01307 01308 return dep_flg; 01309 } 01310 01311 01312 template<class ScalarType, class MV, class OP> 01313 int 01314 DGKSOrthoManager<ScalarType, MV, OP>::blkOrthoSing ( MV &X, Teuchos::RCP<MV> MX, 01315 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 01316 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 01317 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ) const 01318 { 01319 Teuchos::Array<Teuchos::RCP<const MV> > Q (QQ); 01320 01321 const ScalarType ONE = SCT::one(); 01322 const ScalarType ZERO = SCT::zero(); 01323 01324 int nq = Q.size(); 01325 int xc = MVT::GetNumberVecs( X ); 01326 std::vector<int> indX( 1 ); 01327 std::vector<ScalarType> oldDot( 1 ), newDot( 1 ); 01328 01329 std::vector<int> qcs( nq ); 01330 for (int i=0; i<nq; i++) { 01331 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01332 } 01333 01334 // Create pointers for the previous vectors of X that have already been orthonormalized. 01335 Teuchos::RCP<const MV> lastQ; 01336 Teuchos::RCP<MV> Xj, MXj; 01337 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lastC; 01338 01339 // Perform the Gram-Schmidt transformation for each vector in the block of vectors. 01340 for (int j=0; j<xc; j++) { 01341 01342 bool dep_flg = false; 01343 01344 // Get a view of the previously orthogonalized vectors and B, add it to the arrays. 01345 if (j > 0) { 01346 std::vector<int> index( j ); 01347 for (int ind=0; ind<j; ind++) { 01348 index[ind] = ind; 01349 } 01350 lastQ = MVT::CloneView( X, index ); 01351 01352 // Add these views to the Q and C arrays. 01353 Q.push_back( lastQ ); 01354 C.push_back( B ); 01355 qcs.push_back( MVT::GetNumberVecs( *lastQ ) ); 01356 } 01357 01358 // Get a view of the current vector in X to orthogonalize. 01359 indX[0] = j; 01360 Xj = MVT::CloneViewNonConst( X, indX ); 01361 if (this->_hasOp) { 01362 MXj = MVT::CloneViewNonConst( *MX, indX ); 01363 } 01364 else { 01365 MXj = Xj; 01366 } 01367 01368 // Compute the initial Op-norms 01369 MVT::MvDot( *Xj, *MXj, oldDot ); 01370 01371 Teuchos::Array<Teuchos::RCP<MV> > MQ(Q.size()); 01372 // Define the product Q^T * (Op*X) 01373 for (int i=0; i<Q.size(); i++) { 01374 01375 // Get a view of the current serial dense matrix 01376 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j ); 01377 01378 // Multiply Q' with MX 01379 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*Xj,MXj,tempC); 01380 // Multiply by Q and subtract the result in Xj 01381 MVT::MvTimesMatAddMv( -ONE, *Q[i], tempC, ONE, *Xj ); 01382 01383 // Update MXj, with the least number of applications of Op as possible 01384 if (this->_hasOp) { 01385 if (xc <= qcs[i]) { 01386 OPT::Apply( *(this->_Op), *Xj, *MXj); 01387 } 01388 else { 01389 // this will possibly be used again below; don't delete it 01390 MQ[i] = MVT::Clone( *Q[i], qcs[i] ); 01391 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] ); 01392 MVT::MvTimesMatAddMv( -ONE, *MQ[i], tempC, ONE, *MXj ); 01393 } 01394 } 01395 } 01396 01397 // Compute the Op-norms 01398 MVT::MvDot( *Xj, *MXj, newDot ); 01399 01400 // Do one step of classical Gram-Schmidt orthogonalization 01401 // with a second correction step if needed. 01402 01403 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*dep_tol_) ) { 01404 01405 for (int i=0; i<Q.size(); i++) { 01406 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j ); 01407 Teuchos::SerialDenseMatrix<int,ScalarType> C2( qcs[i], 1 ); 01408 01409 // Apply another step of classical Gram-Schmidt 01410 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*Xj,MXj,C2); 01411 tempC += C2; 01412 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, *Xj ); 01413 01414 // Update MXj, with the least number of applications of Op as possible 01415 if (this->_hasOp) { 01416 if (MQ[i].get()) { 01417 // MQ was allocated and computed above; use it 01418 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MXj ); 01419 } 01420 else if (xc <= qcs[i]) { 01421 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01422 OPT::Apply( *(this->_Op), *Xj, *MXj); 01423 } 01424 } 01425 } // for (int i=0; i<Q.size(); i++) 01426 01427 // Compute the Op-norms after the correction step. 01428 MVT::MvDot( *Xj, *MXj, newDot ); 01429 01430 } // if () 01431 01432 // Check for linear dependence. 01433 if (SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*sing_tol_)) { 01434 dep_flg = true; 01435 } 01436 01437 // Normalize the new vector if it's not dependent 01438 if (!dep_flg) { 01439 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01440 01441 MVT::MvAddMv( ONE/diag, *Xj, ZERO, *Xj, *Xj ); 01442 if (this->_hasOp) { 01443 // Update MXj. 01444 MVT::MvAddMv( ONE/diag, *MXj, ZERO, *MXj, *MXj ); 01445 } 01446 01447 // Enter value on diagonal of B. 01448 (*B)(j,j) = diag; 01449 } 01450 else { 01451 // Create a random vector and orthogonalize it against all previous columns of Q. 01452 Teuchos::RCP<MV> tempXj = MVT::Clone( X, 1 ); 01453 Teuchos::RCP<MV> tempMXj; 01454 MVT::MvRandom( *tempXj ); 01455 if (this->_hasOp) { 01456 tempMXj = MVT::Clone( X, 1 ); 01457 OPT::Apply( *(this->_Op), *tempXj, *tempMXj ); 01458 } 01459 else { 01460 tempMXj = tempXj; 01461 } 01462 MVT::MvDot( *tempXj, *tempMXj, oldDot ); 01463 // 01464 for (int num_orth=0; num_orth<max_blk_ortho_; num_orth++) { 01465 01466 for (int i=0; i<Q.size(); i++) { 01467 Teuchos::SerialDenseMatrix<int,ScalarType> product( qcs[i], 1 ); 01468 01469 // Apply another step of classical Gram-Schmidt 01470 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*tempXj,tempMXj,product); 01471 MVT::MvTimesMatAddMv( -ONE, *Q[i], product, ONE, *tempXj ); 01472 01473 // Update MXj, with the least number of applications of Op as possible 01474 if (this->_hasOp) { 01475 if (MQ[i].get()) { 01476 // MQ was allocated and computed above; use it 01477 MVT::MvTimesMatAddMv( -ONE, *MQ[i], product, ONE, *tempMXj ); 01478 } 01479 else if (xc <= qcs[i]) { 01480 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01481 OPT::Apply( *(this->_Op), *tempXj, *tempMXj); 01482 } 01483 } 01484 } // for (int i=0; i<nq; i++) 01485 01486 } 01487 01488 // Compute the Op-norms after the correction step. 01489 MVT::MvDot( *tempXj, *tempMXj, newDot ); 01490 01491 // Copy vector into current column of Xj 01492 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) { 01493 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01494 01495 // Enter value on diagonal of B. 01496 (*B)(j,j) = ZERO; 01497 01498 // Copy vector into current column of _basisvecs 01499 MVT::MvAddMv( ONE/diag, *tempXj, ZERO, *tempXj, *Xj ); 01500 if (this->_hasOp) { 01501 MVT::MvAddMv( ONE/diag, *tempMXj, ZERO, *tempMXj, *MXj ); 01502 } 01503 } 01504 else { 01505 return j; 01506 } 01507 } // if (!dep_flg) 01508 01509 // Remove the vectors from array 01510 if (j > 0) { 01511 Q.resize( nq ); 01512 C.resize( nq ); 01513 qcs.resize( nq ); 01514 } 01515 01516 } // for (int j=0; j<xc; j++) 01517 01518 return xc; 01519 } 01520 01521 } // namespace Belos 01522 01523 #endif // BELOS_DGKS_ORTHOMANAGER_HPP 01524
1.7.6.1