|
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_ICGS_ORTHOMANAGER_HPP 00048 #define BELOS_ICGS_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 getDefaultIcgsParameters() 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 ICGS orthogonalization. 00085 // Documentation will be embedded in the parameter list. 00086 const int defaultMaxNumOrthogPasses = 2; 00087 const magnitude_type eps = STM::eps(); 00088 const magnitude_type defaultBlkTol = 00089 as<magnitude_type> (10) * STM::squareroot (eps); 00090 const magnitude_type defaultSingTol = as<magnitude_type> (10) * eps; 00091 00092 Teuchos::RCP<Teuchos::ParameterList> params = 00093 Teuchos::parameterList ("ICGS"); 00094 params->set ("maxNumOrthogPasses", defaultMaxNumOrthogPasses, 00095 "Maximum number of orthogonalization passes " 00096 "(includes the first). Default is 2, since " 00097 "\"twice is enough\" for Krylov methods."); 00098 params->set ("blkTol", defaultBlkTol, 00099 "Block reorthogonalization threshhold."); 00100 params->set ("singTol", defaultSingTol, 00101 "Singular block detection threshold."); 00102 return params; 00103 } 00104 00109 template<class ScalarType> 00110 Teuchos::RCP<const Teuchos::ParameterList> TEUCHOS_DEPRECATED 00111 getFastIcgsParameters() 00112 { 00113 using Teuchos::ParameterList; 00114 using Teuchos::RCP; 00115 using Teuchos::rcp; 00116 using Teuchos::ScalarTraits; 00117 typedef typename ScalarTraits<ScalarType>::magnitudeType magnitude_type; 00118 typedef ScalarTraits<magnitude_type> STM; 00119 00120 RCP<const ParameterList> defaultParams = getDefaultIcgsParameters<ScalarType>(); 00121 // Start with a clone of the default parameters. 00122 RCP<ParameterList> params = rcp (new ParameterList (*defaultParams)); 00123 00124 const int maxBlkOrtho = 1; 00125 const magnitude_type blkTol = STM::zero(); 00126 const magnitude_type singTol = STM::zero(); 00127 00128 params->set ("maxNumOrthogPasses", maxBlkOrtho); 00129 params->set ("blkTol", blkTol); 00130 params->set ("singTol", singTol); 00131 00132 return params; 00133 } 00134 00144 template<class ScalarType> 00145 void TEUCHOS_DEPRECATED 00146 readIcgsParameters (const Teuchos::RCP<const Teuchos::ParameterList>& params, 00147 int& maxNumOrthogPasses, 00148 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType& blkTol, 00149 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType& singTol) 00150 { 00151 using Teuchos::ParameterList; 00152 using Teuchos::RCP; 00153 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType magnitude_type; 00154 const magnitude_type zero = Teuchos::ScalarTraits<magnitude_type>::zero(); 00155 RCP<const ParameterList> defaultParams = getDefaultIcgsParameters<ScalarType>(); 00156 00157 // Using temporary variables and fetching all values before 00158 // setting the output arguments ensures the strong exception 00159 // guarantee for this function: if an exception is thrown, no 00160 // externally visible side effects (in this case, setting the 00161 // output arguments) have taken place. 00162 int _maxNumOrthogPasses; 00163 magnitude_type _blkTol, _singTol; 00164 if (params.is_null()) { 00165 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00166 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00167 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00168 } else { 00169 try { 00170 _maxNumOrthogPasses = params->get<int> ("maxNumOrthogPasses"); 00171 if (_maxNumOrthogPasses < 1) 00172 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00173 } catch (Teuchos::Exceptions::InvalidParameter&) { 00174 _maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00175 } 00176 00177 try { 00178 _blkTol = params->get<magnitude_type> ("blkTol"); 00179 if (_blkTol < zero) 00180 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00181 } catch (Teuchos::Exceptions::InvalidParameter&) { 00182 try { 00183 // People may have used depTol instead of blkTol for this 00184 // parameter's name, by analogy with DGKS. 00185 _blkTol = params->get<magnitude_type> ("depTol"); 00186 } catch (Teuchos::Exceptions::InvalidParameter&) { 00187 _blkTol = defaultParams->get<magnitude_type> ("blkTol"); 00188 } 00189 } 00190 00191 try { 00192 _singTol = params->get<magnitude_type> ("singTol"); 00193 if (_singTol < zero) 00194 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00195 } catch (Teuchos::Exceptions::InvalidParameter&) { 00196 _singTol = defaultParams->get<magnitude_type> ("singTol"); 00197 } 00198 } 00199 maxNumOrthogPasses = _maxNumOrthogPasses; 00200 blkTol = _blkTol; 00201 singTol = _singTol; 00202 } 00203 00204 00205 template<class ScalarType, class MV, class OP> 00206 class ICGSOrthoManager : 00207 public MatOrthoManager<ScalarType,MV,OP>, 00208 public Teuchos::ParameterListAcceptorDefaultBase 00209 { 00210 private: 00211 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType; 00212 typedef typename Teuchos::ScalarTraits<MagnitudeType> MGT; 00213 typedef Teuchos::ScalarTraits<ScalarType> SCT; 00214 typedef MultiVecTraits<ScalarType,MV> MVT; 00215 typedef MultiVecTraitsExt<ScalarType,MV> MVText; 00216 typedef OperatorTraits<ScalarType,MV,OP> OPT; 00217 00218 public: 00220 00221 00223 ICGSOrthoManager( const std::string& label = "Belos", 00224 Teuchos::RCP<const OP> Op = Teuchos::null, 00225 const int max_ortho_steps = 2, 00226 const MagnitudeType blk_tol = 10*MGT::squareroot( MGT::eps() ), 00227 const MagnitudeType sing_tol = 10*MGT::eps() ) 00228 : MatOrthoManager<ScalarType,MV,OP>(Op), 00229 max_ortho_steps_( max_ortho_steps ), 00230 blk_tol_( blk_tol ), 00231 sing_tol_( sing_tol ), 00232 label_( label ) 00233 { 00234 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00235 std::string orthoLabel = label_ + ": Orthogonalization"; 00236 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel); 00237 00238 std::string updateLabel = label_ + ": Ortho (Update)"; 00239 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel); 00240 00241 std::string normLabel = label_ + ": Ortho (Norm)"; 00242 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel); 00243 00244 std::string ipLabel = label_ + ": Ortho (Inner Product)"; 00245 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel); 00246 #endif 00247 } 00248 00250 ICGSOrthoManager (const Teuchos::RCP<Teuchos::ParameterList>& plist, 00251 const std::string& label = "Belos", 00252 Teuchos::RCP<const OP> Op = Teuchos::null) : 00253 MatOrthoManager<ScalarType,MV,OP>(Op), 00254 max_ortho_steps_ (2), 00255 blk_tol_ (10 * MGT::squareroot (MGT::eps())), 00256 sing_tol_ (10 * MGT::eps()), 00257 label_ (label) 00258 { 00259 setParameterList (plist); 00260 00261 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00262 std::string orthoLabel = label_ + ": Orthogonalization"; 00263 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel); 00264 00265 std::string updateLabel = label_ + ": Ortho (Update)"; 00266 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel); 00267 00268 std::string normLabel = label_ + ": Ortho (Norm)"; 00269 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel); 00270 00271 std::string ipLabel = label_ + ": Ortho (Inner Product)"; 00272 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel); 00273 #endif 00274 } 00275 00277 ~ICGSOrthoManager() {} 00279 00281 00282 00283 void 00284 setParameterList (const Teuchos::RCP<Teuchos::ParameterList>& plist) 00285 { 00286 using Teuchos::Exceptions::InvalidParameterName; 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 params = parameterList (*defaultParams); 00295 } else { 00296 params = plist; 00297 // Some users might want to specify "blkTol" as "depTol". Due 00298 // to this case, we don't invoke 00299 // validateParametersAndSetDefaults on params. Instead, we go 00300 // through the parameter list one parameter at a time and look 00301 // for alternatives. 00302 } 00303 00304 // Using temporary variables and fetching all values before 00305 // setting the output arguments ensures the strong exception 00306 // guarantee for this function: if an exception is thrown, no 00307 // externally visible side effects (in this case, setting the 00308 // output arguments) have taken place. 00309 int maxNumOrthogPasses; 00310 MagnitudeType blkTol; 00311 MagnitudeType singTol; 00312 00313 try { 00314 maxNumOrthogPasses = params->get<int> ("maxNumOrthogPasses"); 00315 } catch (InvalidParameterName&) { 00316 maxNumOrthogPasses = defaultParams->get<int> ("maxNumOrthogPasses"); 00317 params->set ("maxNumOrthogPasses", maxNumOrthogPasses); 00318 } 00319 00320 // Handling of the "blkTol" parameter is a special case. This 00321 // is because some users may prefer to call this parameter 00322 // "depTol" for consistency with DGKS. However, our default 00323 // parameter list calls this "blkTol", and we don't want the 00324 // default list's value to override the user's value. Thus, we 00325 // first check the user's parameter list for both names, and 00326 // only then access the default parameter list. 00327 try { 00328 blkTol = params->get<MagnitudeType> ("blkTol"); 00329 } catch (InvalidParameterName&) { 00330 try { 00331 blkTol = params->get<MagnitudeType> ("depTol"); 00332 // "depTol" is the wrong name, so remove it and replace with 00333 // "blkTol". We'll set "blkTol" below. 00334 params->remove ("depTol"); 00335 } catch (InvalidParameterName&) { 00336 blkTol = defaultParams->get<MagnitudeType> ("blkTol"); 00337 } 00338 params->set ("blkTol", blkTol); 00339 } 00340 00341 try { 00342 singTol = params->get<MagnitudeType> ("singTol"); 00343 } catch (InvalidParameterName&) { 00344 singTol = defaultParams->get<MagnitudeType> ("singTol"); 00345 params->set ("singTol", singTol); 00346 } 00347 00348 max_ortho_steps_ = maxNumOrthogPasses; 00349 blk_tol_ = blkTol; 00350 sing_tol_ = singTol; 00351 00352 setMyParamList (params); 00353 } 00354 00355 Teuchos::RCP<const Teuchos::ParameterList> 00356 getValidParameters () const 00357 { 00358 using Teuchos::as; 00359 using Teuchos::ParameterList; 00360 using Teuchos::parameterList; 00361 using Teuchos::RCP; 00362 00363 if (defaultParams_.is_null()) { 00364 RCP<ParameterList> params = parameterList ("ICGS"); 00365 00366 // Default parameter values for ICGS orthogonalization. 00367 // Documentation will be embedded in the parameter list. 00368 const int defaultMaxNumOrthogPasses = 2; 00369 const MagnitudeType eps = MGT::eps(); 00370 const MagnitudeType defaultBlkTol = 00371 as<MagnitudeType> (10) * MGT::squareroot (eps); 00372 const MagnitudeType defaultSingTol = as<MagnitudeType> (10) * eps; 00373 00374 params->set ("maxNumOrthogPasses", defaultMaxNumOrthogPasses, 00375 "Maximum number of orthogonalization passes (includes the " 00376 "first). Default is 2, since \"twice is enough\" for Krylov " 00377 "methods."); 00378 params->set ("blkTol", defaultBlkTol, "Block reorthogonalization " 00379 "threshhold."); 00380 params->set ("singTol", defaultSingTol, "Singular block detection " 00381 "threshold."); 00382 defaultParams_ = params; 00383 } 00384 return defaultParams_; 00385 } 00387 00392 Teuchos::RCP<const Teuchos::ParameterList> 00393 getFastParameters () const 00394 { 00395 using Teuchos::as; 00396 using Teuchos::ParameterList; 00397 using Teuchos::parameterList; 00398 using Teuchos::RCP; 00399 00400 RCP<const ParameterList> defaultParams = getValidParameters (); 00401 // Start with a clone of the default parameters. 00402 RCP<ParameterList> params = parameterList (*defaultParams); 00403 00404 const int maxBlkOrtho = 1; // No block reorthogonalization 00405 const MagnitudeType blkTol = MGT::zero(); 00406 const MagnitudeType singTol = MGT::zero(); 00407 00408 params->set ("maxNumOrthogPasses", maxBlkOrtho); 00409 params->set ("blkTol", blkTol); 00410 params->set ("singTol", singTol); 00411 00412 return params; 00413 } 00414 00416 00417 00419 void setBlkTol( const MagnitudeType blk_tol ) { 00420 // Update the parameter list as well. 00421 Teuchos::RCP<Teuchos::ParameterList> params = getNonconstParameterList(); 00422 if (! params.is_null()) { 00423 // If it's null, then we haven't called setParameterList() 00424 // yet. It's entirely possible to construct the parameter 00425 // list on demand, so we don't try to create the parameter 00426 // list here. 00427 params->set ("blkTol", blk_tol); 00428 } 00429 blk_tol_ = blk_tol; 00430 } 00431 00433 void setSingTol( const MagnitudeType sing_tol ) { 00434 // Update the parameter list as well. 00435 Teuchos::RCP<Teuchos::ParameterList> params = getNonconstParameterList(); 00436 if (! params.is_null()) { 00437 // If it's null, then we haven't called setParameterList() 00438 // yet. It's entirely possible to construct the parameter 00439 // list on demand, so we don't try to create the parameter 00440 // list here. 00441 params->set ("singTol", sing_tol); 00442 } 00443 sing_tol_ = sing_tol; 00444 } 00445 00447 MagnitudeType getBlkTol() const { return blk_tol_; } 00448 00450 MagnitudeType getSingTol() const { return sing_tol_; } 00451 00453 00454 00456 00457 00485 void project ( MV &X, Teuchos::RCP<MV> MX, 00486 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00487 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00488 00489 00492 void project ( MV &X, 00493 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00494 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const { 00495 project(X,Teuchos::null,C,Q); 00496 } 00497 00498 00499 00524 int normalize ( MV &X, Teuchos::RCP<MV> MX, 00525 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B) const; 00526 00527 00530 int normalize ( MV &X, Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const { 00531 return normalize(X,Teuchos::null,B); 00532 } 00533 00534 protected: 00535 00577 virtual int 00578 projectAndNormalizeWithMxImpl (MV &X, 00579 Teuchos::RCP<MV> MX, 00580 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00581 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00582 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00583 00584 public: 00586 00587 00588 00592 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00593 orthonormError(const MV &X) const { 00594 return orthonormError(X,Teuchos::null); 00595 } 00596 00601 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00602 orthonormError(const MV &X, Teuchos::RCP<const MV> MX) const; 00603 00607 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00608 orthogError(const MV &X1, const MV &X2) const { 00609 return orthogError(X1,Teuchos::null,X2); 00610 } 00611 00616 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00617 orthogError(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const; 00618 00620 00622 00623 00626 void setLabel(const std::string& label); 00627 00630 const std::string& getLabel() const { return label_; } 00631 00633 00634 private: 00636 int max_ortho_steps_; 00638 MagnitudeType blk_tol_; 00640 MagnitudeType sing_tol_; 00641 00643 std::string label_; 00644 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00645 Teuchos::RCP<Teuchos::Time> timerOrtho_, timerUpdate_, 00646 timerNorm_, timerScale_, timerInnerProd_; 00647 #endif // BELOS_TEUCHOS_TIME_MONITOR 00648 00650 mutable Teuchos::RCP<Teuchos::ParameterList> defaultParams_; 00651 00653 int findBasis(MV &X, Teuchos::RCP<MV> MX, 00654 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > C, 00655 bool completeBasis, int howMany = -1 ) const; 00656 00658 bool blkOrtho1 ( MV &X, Teuchos::RCP<MV> MX, 00659 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00660 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00661 00663 bool blkOrtho ( MV &X, Teuchos::RCP<MV> MX, 00664 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00665 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const; 00666 00680 int blkOrthoSing ( MV &X, Teuchos::RCP<MV> MX, 00681 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00682 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00683 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ) const; 00684 }; 00685 00687 // Set the label for this orthogonalization manager and create new timers if it's changed 00688 template<class ScalarType, class MV, class OP> 00689 void ICGSOrthoManager<ScalarType,MV,OP>::setLabel(const std::string& label) 00690 { 00691 if (label != label_) { 00692 label_ = label; 00693 std::string orthoLabel = label_ + ": Orthogonalization"; 00694 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00695 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel); 00696 #endif 00697 00698 std::string updateLabel = label_ + ": Ortho (Update)"; 00699 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00700 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel); 00701 #endif 00702 00703 std::string normLabel = label_ + ": Ortho (Norm)"; 00704 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00705 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel); 00706 #endif 00707 00708 std::string ipLabel = label_ + ": Ortho (Inner Product)"; 00709 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00710 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel); 00711 #endif 00712 } 00713 } 00714 00716 // Compute the distance from orthonormality 00717 template<class ScalarType, class MV, class OP> 00718 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00719 ICGSOrthoManager<ScalarType,MV,OP>::orthonormError(const MV &X, Teuchos::RCP<const MV> MX) const { 00720 const ScalarType ONE = SCT::one(); 00721 int rank = MVT::GetNumberVecs(X); 00722 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank); 00723 MatOrthoManager<ScalarType,MV,OP>::innerProd(X,X,MX,xTx); 00724 for (int i=0; i<rank; i++) { 00725 xTx(i,i) -= ONE; 00726 } 00727 return xTx.normFrobenius(); 00728 } 00729 00731 // Compute the distance from orthogonality 00732 template<class ScalarType, class MV, class OP> 00733 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType 00734 ICGSOrthoManager<ScalarType,MV,OP>::orthogError(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const { 00735 int r1 = MVT::GetNumberVecs(X1); 00736 int r2 = MVT::GetNumberVecs(X2); 00737 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r2,r1); 00738 MatOrthoManager<ScalarType,MV,OP>::innerProd(X2,X1,MX1,xTx); 00739 return xTx.normFrobenius(); 00740 } 00741 00743 // Find an Op-orthonormal basis for span(X) - span(W) 00744 template<class ScalarType, class MV, class OP> 00745 int 00746 ICGSOrthoManager<ScalarType, MV, OP>:: 00747 projectAndNormalizeWithMxImpl (MV &X, 00748 Teuchos::RCP<MV> MX, 00749 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00750 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 00751 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const 00752 { 00753 using Teuchos::Array; 00754 using Teuchos::null; 00755 using Teuchos::is_null; 00756 using Teuchos::RCP; 00757 using Teuchos::rcp; 00758 using Teuchos::SerialDenseMatrix; 00759 typedef SerialDenseMatrix< int, ScalarType > serial_dense_matrix_type; 00760 typedef typename Array< RCP< const MV > >::size_type size_type; 00761 00762 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00763 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00764 #endif 00765 00766 ScalarType ONE = SCT::one(); 00767 //ScalarType ZERO = SCT::zero(); 00768 00769 int nq = Q.size(); 00770 int xc = MVT::GetNumberVecs( X ); 00771 ptrdiff_t xr = MVText::GetGlobalLength( X ); 00772 int rank = xc; 00773 00774 // If the user doesn't want to store the normalization 00775 // coefficients, allocate some local memory for them. This will 00776 // go away at the end of this method. 00777 if (is_null (B)) { 00778 B = rcp (new serial_dense_matrix_type (xc, xc)); 00779 } 00780 // Likewise, if the user doesn't want to store the projection 00781 // coefficients, allocate some local memory for them. Also make 00782 // sure that all the entries of C are the right size. We're going 00783 // to overwrite them anyway, so we don't have to worry about the 00784 // contents (other than to resize them if they are the wrong 00785 // size). 00786 if (C.size() < nq) 00787 C.resize (nq); 00788 for (size_type k = 0; k < nq; ++k) 00789 { 00790 const int numRows = MVT::GetNumberVecs (*Q[k]); 00791 const int numCols = xc; // Number of vectors in X 00792 00793 if (is_null (C[k])) 00794 C[k] = rcp (new serial_dense_matrix_type (numRows, numCols)); 00795 else if (C[k]->numRows() != numRows || C[k]->numCols() != numCols) 00796 { 00797 int err = C[k]->reshape (numRows, numCols); 00798 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, 00799 "IMGS orthogonalization: failed to reshape " 00800 "C[" << k << "] (the array of block " 00801 "coefficients resulting from projecting X " 00802 "against Q[1:" << nq << "])."); 00803 } 00804 } 00805 00806 /****** DO NOT MODIFY *MX IF _hasOp == false ******/ 00807 if (this->_hasOp) { 00808 if (MX == Teuchos::null) { 00809 // we need to allocate space for MX 00810 MX = MVT::Clone(X,MVT::GetNumberVecs(X)); 00811 OPT::Apply(*(this->_Op),X,*MX); 00812 } 00813 } 00814 else { 00815 // Op == I --> MX = X (ignore it if the user passed it in) 00816 MX = Teuchos::rcp( &X, false ); 00817 } 00818 00819 int mxc = MVT::GetNumberVecs( *MX ); 00820 ptrdiff_t mxr = MVText::GetGlobalLength( *MX ); 00821 00822 // short-circuit 00823 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument, "Belos::ICGSOrthoManager::projectAndNormalize(): X must be non-empty" ); 00824 00825 int numbas = 0; 00826 for (int i=0; i<nq; i++) { 00827 numbas += MVT::GetNumberVecs( *Q[i] ); 00828 } 00829 00830 // check size of B 00831 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument, 00832 "Belos::ICGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of B" ); 00833 // check size of X and MX 00834 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument, 00835 "Belos::ICGSOrthoManager::projectAndNormalize(): MVT returned negative dimensions for X,MX" ); 00836 // check size of X w.r.t. MX 00837 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument, 00838 "Belos::ICGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of MX" ); 00839 // check feasibility 00840 //TEUCHOS_TEST_FOR_EXCEPTION( numbas+xc > xr, std::invalid_argument, 00841 // "Belos::ICGSOrthoManager::projectAndNormalize(): Orthogonality constraints not feasible" ); 00842 00843 // Some flags for checking dependency returns from the internal orthogonalization methods 00844 bool dep_flg = false; 00845 00846 // Make a temporary copy of X and MX, just in case a block dependency is detected. 00847 Teuchos::RCP<MV> tmpX, tmpMX; 00848 tmpX = MVT::CloneCopy(X); 00849 if (this->_hasOp) { 00850 tmpMX = MVT::CloneCopy(*MX); 00851 } 00852 00853 if (xc == 1) { 00854 00855 // Use the cheaper block orthogonalization. 00856 // NOTE: Don't check for dependencies because the update has one vector. 00857 dep_flg = blkOrtho1( X, MX, C, Q ); 00858 00859 // Normalize the new block X 00860 { 00861 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00862 Teuchos::TimeMonitor normTimer( *timerNorm_ ); 00863 #endif 00864 if ( B == Teuchos::null ) { 00865 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) ); 00866 } 00867 std::vector<ScalarType> diag(xc); 00868 MVT::MvDot( X, *MX, diag ); 00869 (*B)(0,0) = SCT::squareroot(SCT::magnitude(diag[0])); 00870 rank = 1; 00871 MVT::MvScale( X, ONE/(*B)(0,0) ); 00872 if (this->_hasOp) { 00873 // Update MXj. 00874 MVT::MvScale( *MX, ONE/(*B)(0,0) ); 00875 } 00876 } 00877 00878 } 00879 else { 00880 00881 // Use the cheaper block orthogonalization. 00882 dep_flg = blkOrtho( X, MX, C, Q ); 00883 00884 // If a dependency has been detected in this block, then perform 00885 // the more expensive nonblock (single vector at a time) 00886 // orthogonalization. 00887 if (dep_flg) { 00888 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q ); 00889 00890 // Copy tmpX back into X. 00891 MVT::Assign( *tmpX, X ); 00892 if (this->_hasOp) { 00893 MVT::Assign( *tmpMX, *MX ); 00894 } 00895 } 00896 else { 00897 // There is no dependency, so orthonormalize new block X 00898 rank = findBasis( X, MX, B, false ); 00899 if (rank < xc) { 00900 // A dependency was found during orthonormalization of X, 00901 // rerun orthogonalization using more expensive nonblock 00902 // orthogonalization. 00903 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q ); 00904 00905 // Copy tmpX back into X. 00906 MVT::Assign( *tmpX, X ); 00907 if (this->_hasOp) { 00908 MVT::Assign( *tmpMX, *MX ); 00909 } 00910 } 00911 } 00912 } // if (xc == 1) { 00913 00914 // this should not raise an std::exception; but our post-conditions oblige us to check 00915 TEUCHOS_TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error, 00916 "Belos::ICGSOrthoManager::projectAndNormalize(): Debug error in rank variable." ); 00917 00918 // Return the rank of X. 00919 return rank; 00920 } 00921 00922 00923 00925 // Find an Op-orthonormal basis for span(X), with rank numvectors(X) 00926 template<class ScalarType, class MV, class OP> 00927 int ICGSOrthoManager<ScalarType, MV, OP>::normalize( 00928 MV &X, Teuchos::RCP<MV> MX, 00929 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const { 00930 00931 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00932 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00933 #endif 00934 00935 // call findBasis, with the instruction to try to generate a basis of rank numvecs(X) 00936 return findBasis(X, MX, B, true); 00937 } 00938 00939 00940 00942 template<class ScalarType, class MV, class OP> 00943 void ICGSOrthoManager<ScalarType, MV, OP>::project( 00944 MV &X, Teuchos::RCP<MV> MX, 00945 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 00946 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const { 00947 // For the inner product defined by the operator Op or the identity (Op == 0) 00948 // -> Orthogonalize X against each Q[i] 00949 // Modify MX accordingly 00950 // 00951 // Note that when Op is 0, MX is not referenced 00952 // 00953 // Parameter variables 00954 // 00955 // X : Vectors to be transformed 00956 // 00957 // MX : Image of the block of vectors X by the mass matrix 00958 // 00959 // Q : Bases to orthogonalize against. These are assumed orthonormal, mutually and independently. 00960 // 00961 00962 #ifdef BELOS_TEUCHOS_TIME_MONITOR 00963 Teuchos::TimeMonitor orthotimer(*timerOrtho_); 00964 #endif 00965 00966 int xc = MVT::GetNumberVecs( X ); 00967 ptrdiff_t xr = MVText::GetGlobalLength( X ); 00968 int nq = Q.size(); 00969 std::vector<int> qcs(nq); 00970 // short-circuit 00971 if (nq == 0 || xc == 0 || xr == 0) { 00972 return; 00973 } 00974 ptrdiff_t qr = MVText::GetGlobalLength ( *Q[0] ); 00975 // if we don't have enough C, expand it with null references 00976 // if we have too many, resize to throw away the latter ones 00977 // if we have exactly as many as we have Q, this call has no effect 00978 C.resize(nq); 00979 00980 00981 /****** DO NOT MODIFY *MX IF _hasOp == false ******/ 00982 if (this->_hasOp) { 00983 if (MX == Teuchos::null) { 00984 // we need to allocate space for MX 00985 MX = MVT::Clone(X,MVT::GetNumberVecs(X)); 00986 OPT::Apply(*(this->_Op),X,*MX); 00987 } 00988 } 00989 else { 00990 // Op == I --> MX = X (ignore it if the user passed it in) 00991 MX = Teuchos::rcp( &X, false ); 00992 } 00993 int mxc = MVT::GetNumberVecs( *MX ); 00994 ptrdiff_t mxr = MVText::GetGlobalLength( *MX ); 00995 00996 // check size of X and Q w.r.t. common sense 00997 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument, 00998 "Belos::ICGSOrthoManager::project(): MVT returned negative dimensions for X,MX" ); 00999 // check size of X w.r.t. MX and Q 01000 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr || xr!=qr, std::invalid_argument, 01001 "Belos::ICGSOrthoManager::project(): Size of X not consistant with MX,Q" ); 01002 01003 // tally up size of all Q and check/allocate C 01004 int baslen = 0; 01005 for (int i=0; i<nq; i++) { 01006 TEUCHOS_TEST_FOR_EXCEPTION( MVText::GetGlobalLength( *Q[i] ) != qr, std::invalid_argument, 01007 "Belos::ICGSOrthoManager::project(): Q lengths not mutually consistant" ); 01008 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01009 TEUCHOS_TEST_FOR_EXCEPTION( qr < qcs[i], std::invalid_argument, 01010 "Belos::ICGSOrthoManager::project(): Q has less rows than columns" ); 01011 baslen += qcs[i]; 01012 01013 // check size of C[i] 01014 if ( C[i] == Teuchos::null ) { 01015 C[i] = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) ); 01016 } 01017 else { 01018 TEUCHOS_TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc , std::invalid_argument, 01019 "Belos::ICGSOrthoManager::project(): Size of Q not consistant with size of C" ); 01020 } 01021 } 01022 01023 // Use the cheaper block orthogonalization, don't check for rank deficiency. 01024 blkOrtho( X, MX, C, Q ); 01025 01026 } 01027 01029 // Find an Op-orthonormal basis for span(X), with the option of extending the subspace so that 01030 // the rank is numvectors(X) 01031 template<class ScalarType, class MV, class OP> 01032 int 01033 ICGSOrthoManager<ScalarType, MV, OP>:: 01034 findBasis (MV &X, 01035 Teuchos::RCP<MV> MX, 01036 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 01037 bool completeBasis, 01038 int howMany) const 01039 { 01040 // For the inner product defined by the operator Op or the identity (Op == 0) 01041 // -> Orthonormalize X 01042 // Modify MX accordingly 01043 // 01044 // Note that when Op is 0, MX is not referenced 01045 // 01046 // Parameter variables 01047 // 01048 // X : Vectors to be orthonormalized 01049 // MX : Image of the multivector X under the operator Op 01050 // Op : Pointer to the operator for the inner product 01051 // 01052 using Teuchos::as; 01053 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01054 Teuchos::TimeMonitor normTimer (*timerNorm_); 01055 #endif // BELOS_TEUCHOS_TIME_MONITOR 01056 01057 const ScalarType ONE = SCT::one (); 01058 const MagnitudeType ZERO = SCT::magnitude (SCT::zero ()); 01059 01060 const int xc = MVT::GetNumberVecs (X); 01061 const ptrdiff_t xr = MVText::GetGlobalLength (X); 01062 01063 if (howMany == -1) { 01064 howMany = xc; 01065 } 01066 01067 /******************************************************* 01068 * If _hasOp == false, we will not reference MX below * 01069 *******************************************************/ 01070 01071 // if Op==null, MX == X (via pointer) 01072 // Otherwise, either the user passed in MX or we will allocated and compute it 01073 if (this->_hasOp) { 01074 if (MX == Teuchos::null) { 01075 // we need to allocate space for MX 01076 MX = MVT::Clone(X,xc); 01077 OPT::Apply(*(this->_Op),X,*MX); 01078 } 01079 } 01080 01081 /* if the user doesn't want to store the coefficienets, 01082 * allocate some local memory for them 01083 */ 01084 if ( B == Teuchos::null ) { 01085 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) ); 01086 } 01087 01088 const int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc; 01089 const ptrdiff_t mxr = (this->_hasOp) ? MVText::GetGlobalLength( *MX ) : xr; 01090 01091 // check size of C, B 01092 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument, 01093 "Belos::ICGSOrthoManager::findBasis(): X must be non-empty" ); 01094 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument, 01095 "Belos::ICGSOrthoManager::findBasis(): Size of X not consistant with size of B" ); 01096 TEUCHOS_TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument, 01097 "Belos::ICGSOrthoManager::findBasis(): Size of X not consistant with size of MX" ); 01098 TEUCHOS_TEST_FOR_EXCEPTION( as<ptrdiff_t> (xc) > xr, std::invalid_argument, 01099 "Belos::ICGSOrthoManager::findBasis(): Size of X not feasible for normalization" ); 01100 TEUCHOS_TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::invalid_argument, 01101 "Belos::ICGSOrthoManager::findBasis(): Invalid howMany parameter" ); 01102 01103 /* xstart is which column we are starting the process with, based on howMany 01104 * columns before xstart are assumed to be Op-orthonormal already 01105 */ 01106 int xstart = xc - howMany; 01107 01108 for (int j = xstart; j < xc; j++) { 01109 01110 // numX is 01111 // * number of currently orthonormal columns of X 01112 // * the index of the current column of X 01113 int numX = j; 01114 bool addVec = false; 01115 01116 // Get a view of the vector currently being worked on. 01117 std::vector<int> index(1); 01118 index[0] = numX; 01119 Teuchos::RCP<MV> Xj = MVT::CloneViewNonConst( X, index ); 01120 Teuchos::RCP<MV> MXj; 01121 if (this->_hasOp) { 01122 // MXj is a view of the current vector in MX 01123 MXj = MVT::CloneViewNonConst( *MX, index ); 01124 } 01125 else { 01126 // MXj is a pointer to Xj, and MUST NOT be modified 01127 MXj = Xj; 01128 } 01129 01130 // Get a view of the previous vectors. 01131 std::vector<int> prev_idx( numX ); 01132 Teuchos::RCP<const MV> prevX, prevMX; 01133 01134 if (numX > 0) { 01135 for (int i=0; i<numX; i++) { 01136 prev_idx[i] = i; 01137 } 01138 prevX = MVT::CloneView( X, prev_idx ); 01139 if (this->_hasOp) { 01140 prevMX = MVT::CloneView( *MX, prev_idx ); 01141 } 01142 } 01143 01144 // Make storage for these Gram-Schmidt iterations. 01145 Teuchos::SerialDenseMatrix<int,ScalarType> product(numX, 1); 01146 std::vector<ScalarType> oldDot( 1 ), newDot( 1 ); 01147 // 01148 // Save old MXj vector and compute Op-norm 01149 // 01150 Teuchos::RCP<MV> oldMXj = MVT::CloneCopy( *MXj ); 01151 MVT::MvDot( *Xj, *MXj, oldDot ); 01152 // Xj^H Op Xj should be real and positive, by the hermitian positive definiteness of Op 01153 TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(oldDot[0]) < ZERO, OrthoError, 01154 "Belos::ICGSOrthoManager::findBasis(): Negative definiteness discovered in inner product" ); 01155 01156 if (numX > 0) { 01157 01158 Teuchos::SerialDenseMatrix<int,ScalarType> P2(numX,1); 01159 01160 for (int i=0; i<max_ortho_steps_; ++i) { 01161 01162 // product <- prevX^T MXj 01163 { 01164 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01165 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01166 #endif 01167 MatOrthoManager<ScalarType,MV,OP>::innerProd(*prevX,*Xj,MXj,P2); 01168 } 01169 01170 // Xj <- Xj - prevX prevX^T MXj 01171 // = Xj - prevX product 01172 { 01173 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01174 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01175 #endif 01176 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj ); 01177 } 01178 01179 // Update MXj 01180 if (this->_hasOp) { 01181 // MXj <- Op*Xj_new 01182 // = Op*(Xj_old - prevX prevX^T MXj) 01183 // = MXj - prevMX product 01184 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01185 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01186 #endif 01187 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj ); 01188 } 01189 01190 // Set coefficients 01191 if ( i==0 ) 01192 product = P2; 01193 else 01194 product += P2; 01195 } 01196 01197 } // if (numX > 0) 01198 01199 // Compute Op-norm with old MXj 01200 MVT::MvDot( *Xj, *oldMXj, newDot ); 01201 01202 // Check to see if the new vector is dependent. 01203 if (completeBasis) { 01204 // 01205 // We need a complete basis, so add random vectors if necessary 01206 // 01207 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(sing_tol_*oldDot[0]) ) { 01208 01209 // Add a random vector and orthogonalize it against previous vectors in block. 01210 addVec = true; 01211 #ifdef ORTHO_DEBUG 01212 std::cout << "Belos::ICGSOrthoManager::findBasis() --> Random for column " << numX << std::endl; 01213 #endif 01214 // 01215 Teuchos::RCP<MV> tempXj = MVT::Clone( X, 1 ); 01216 Teuchos::RCP<MV> tempMXj; 01217 MVT::MvRandom( *tempXj ); 01218 if (this->_hasOp) { 01219 tempMXj = MVT::Clone( X, 1 ); 01220 OPT::Apply( *(this->_Op), *tempXj, *tempMXj ); 01221 } 01222 else { 01223 tempMXj = tempXj; 01224 } 01225 MVT::MvDot( *tempXj, *tempMXj, oldDot ); 01226 // 01227 for (int num_orth=0; num_orth<max_ortho_steps_; num_orth++){ 01228 { 01229 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01230 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01231 #endif 01232 MatOrthoManager<ScalarType,MV,OP>::innerProd(*prevX,*tempXj,tempMXj,product); 01233 } 01234 { 01235 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01236 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01237 #endif 01238 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *tempXj ); 01239 } 01240 if (this->_hasOp) { 01241 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01242 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01243 #endif 01244 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *tempMXj ); 01245 } 01246 } 01247 // Compute new Op-norm 01248 MVT::MvDot( *tempXj, *tempMXj, newDot ); 01249 // 01250 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) { 01251 // Copy vector into current column of _basisvecs 01252 MVT::Assign( *tempXj, *Xj ); 01253 if (this->_hasOp) { 01254 MVT::Assign( *tempMXj, *MXj ); 01255 } 01256 } 01257 else { 01258 return numX; 01259 } 01260 } 01261 } 01262 else { 01263 // 01264 // We only need to detect dependencies. 01265 // 01266 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*blk_tol_) ) { 01267 return numX; 01268 } 01269 } 01270 01271 // If we haven't left this method yet, then we can normalize the new vector Xj. 01272 // Normalize Xj. 01273 // Xj <- Xj / std::sqrt(newDot) 01274 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01275 { 01276 MVT::MvScale( *Xj, ONE/diag ); 01277 if (this->_hasOp) { 01278 // Update MXj. 01279 MVT::MvScale( *MXj, ONE/diag ); 01280 } 01281 } 01282 01283 // If we've added a random vector, enter a zero in the j'th diagonal element. 01284 if (addVec) { 01285 (*B)(j,j) = ZERO; 01286 } 01287 else { 01288 (*B)(j,j) = diag; 01289 } 01290 01291 // Save the coefficients, if we are working on the original vector and not a randomly generated one 01292 if (!addVec) { 01293 for (int i=0; i<numX; i++) { 01294 (*B)(i,j) = product(i,0); 01295 } 01296 } 01297 01298 } // for (j = 0; j < xc; ++j) 01299 01300 return xc; 01301 } 01302 01304 // Routine to compute the block orthogonalization 01305 template<class ScalarType, class MV, class OP> 01306 bool 01307 ICGSOrthoManager<ScalarType, MV, OP>::blkOrtho1 ( MV &X, Teuchos::RCP<MV> MX, 01308 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 01309 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const 01310 { 01311 int nq = Q.size(); 01312 int xc = MVT::GetNumberVecs( X ); 01313 const ScalarType ONE = SCT::one(); 01314 01315 std::vector<int> qcs( nq ); 01316 for (int i=0; i<nq; i++) { 01317 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01318 } 01319 01320 // Perform the Gram-Schmidt transformation for a block of vectors 01321 01322 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq); 01323 // Define the product Q^T * (Op*X) 01324 for (int i=0; i<nq; i++) { 01325 // Multiply Q' with MX 01326 { 01327 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01328 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01329 #endif 01330 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,*C[i]); 01331 } 01332 // Multiply by Q and subtract the result in X 01333 { 01334 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01335 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01336 #endif 01337 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X ); 01338 } 01339 01340 // Update MX, with the least number of applications of Op as possible 01341 if (this->_hasOp) { 01342 if (xc <= qcs[i]) { 01343 OPT::Apply( *(this->_Op), X, *MX); 01344 } 01345 else { 01346 // this will possibly be used again below; don't delete it 01347 MQ[i] = MVT::Clone( *Q[i], qcs[i] ); 01348 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] ); 01349 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX ); 01350 } 01351 } 01352 } 01353 01354 // Do as many steps of classical Gram-Schmidt as required by max_ortho_steps_ 01355 for (int j = 1; j < max_ortho_steps_; ++j) { 01356 01357 for (int i=0; i<nq; i++) { 01358 Teuchos::SerialDenseMatrix<int,ScalarType> C2(*C[i]); 01359 01360 // Apply another step of classical Gram-Schmidt 01361 { 01362 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01363 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01364 #endif 01365 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,C2); 01366 } 01367 *C[i] += C2; 01368 { 01369 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01370 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01371 #endif 01372 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X ); 01373 } 01374 01375 // Update MX, with the least number of applications of Op as possible 01376 if (this->_hasOp) { 01377 if (MQ[i].get()) { 01378 // MQ was allocated and computed above; use it 01379 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX ); 01380 } 01381 else if (xc <= qcs[i]) { 01382 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01383 OPT::Apply( *(this->_Op), X, *MX); 01384 } 01385 } 01386 } // for (int i=0; i<nq; i++) 01387 } // for (int j = 0; j < max_ortho_steps; ++j) 01388 01389 return false; 01390 } 01391 01393 // Routine to compute the block orthogonalization 01394 template<class ScalarType, class MV, class OP> 01395 bool 01396 ICGSOrthoManager<ScalarType, MV, OP>::blkOrtho ( MV &X, Teuchos::RCP<MV> MX, 01397 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 01398 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const 01399 { 01400 int nq = Q.size(); 01401 int xc = MVT::GetNumberVecs( X ); 01402 bool dep_flg = false; 01403 const ScalarType ONE = SCT::one(); 01404 01405 std::vector<int> qcs( nq ); 01406 for (int i=0; i<nq; i++) { 01407 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01408 } 01409 01410 // Perform the Gram-Schmidt transformation for a block of vectors 01411 01412 // Compute the initial Op-norms 01413 std::vector<ScalarType> oldDot( xc ); 01414 MVT::MvDot( X, *MX, oldDot ); 01415 01416 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq); 01417 // Define the product Q^T * (Op*X) 01418 for (int i=0; i<nq; i++) { 01419 // Multiply Q' with MX 01420 { 01421 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01422 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01423 #endif 01424 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,*C[i]); 01425 } 01426 // Multiply by Q and subtract the result in X 01427 { 01428 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01429 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01430 #endif 01431 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X ); 01432 } 01433 // Update MX, with the least number of applications of Op as possible 01434 if (this->_hasOp) { 01435 if (xc <= qcs[i]) { 01436 OPT::Apply( *(this->_Op), X, *MX); 01437 } 01438 else { 01439 // this will possibly be used again below; don't delete it 01440 MQ[i] = MVT::Clone( *Q[i], qcs[i] ); 01441 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] ); 01442 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX ); 01443 } 01444 } 01445 } 01446 01447 // Do as many steps of classical Gram-Schmidt as required by max_ortho_steps_ 01448 for (int j = 1; j < max_ortho_steps_; ++j) { 01449 01450 for (int i=0; i<nq; i++) { 01451 Teuchos::SerialDenseMatrix<int,ScalarType> C2(*C[i]); 01452 01453 // Apply another step of classical Gram-Schmidt 01454 { 01455 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01456 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ ); 01457 #endif 01458 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],X,MX,C2); 01459 } 01460 *C[i] += C2; 01461 { 01462 #ifdef BELOS_TEUCHOS_TIME_MONITOR 01463 Teuchos::TimeMonitor updateTimer( *timerUpdate_ ); 01464 #endif 01465 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X ); 01466 } 01467 01468 // Update MX, with the least number of applications of Op as possible 01469 if (this->_hasOp) { 01470 if (MQ[i].get()) { 01471 // MQ was allocated and computed above; use it 01472 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX ); 01473 } 01474 else if (xc <= qcs[i]) { 01475 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01476 OPT::Apply( *(this->_Op), X, *MX); 01477 } 01478 } 01479 } // for (int i=0; i<nq; i++) 01480 } // for (int j = 0; j < max_ortho_steps; ++j) 01481 01482 // Compute new Op-norms 01483 std::vector<ScalarType> newDot(xc); 01484 MVT::MvDot( X, *MX, newDot ); 01485 01486 // Check to make sure the new block of vectors are not dependent on previous vectors 01487 for (int i=0; i<xc; i++){ 01488 if (SCT::magnitude(newDot[i]) < SCT::magnitude(oldDot[i] * blk_tol_)) { 01489 dep_flg = true; 01490 break; 01491 } 01492 } // end for (i=0;...) 01493 01494 return dep_flg; 01495 } 01496 01497 template<class ScalarType, class MV, class OP> 01498 int 01499 ICGSOrthoManager<ScalarType, MV, OP>::blkOrthoSing ( MV &X, Teuchos::RCP<MV> MX, 01500 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C, 01501 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B, 01502 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ) const 01503 { 01504 Teuchos::Array<Teuchos::RCP<const MV> > Q (QQ); 01505 01506 const ScalarType ONE = SCT::one(); 01507 const ScalarType ZERO = SCT::zero(); 01508 01509 int nq = Q.size(); 01510 int xc = MVT::GetNumberVecs( X ); 01511 std::vector<int> indX( 1 ); 01512 std::vector<ScalarType> oldDot( 1 ), newDot( 1 ); 01513 01514 std::vector<int> qcs( nq ); 01515 for (int i=0; i<nq; i++) { 01516 qcs[i] = MVT::GetNumberVecs( *Q[i] ); 01517 } 01518 01519 // Create pointers for the previous vectors of X that have already been orthonormalized. 01520 Teuchos::RCP<const MV> lastQ; 01521 Teuchos::RCP<MV> Xj, MXj; 01522 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lastC; 01523 01524 // Perform the Gram-Schmidt transformation for each vector in the block of vectors. 01525 for (int j=0; j<xc; j++) { 01526 01527 bool dep_flg = false; 01528 01529 // Get a view of the previously orthogonalized vectors and B, add it to the arrays. 01530 if (j > 0) { 01531 std::vector<int> index( j ); 01532 for (int ind=0; ind<j; ind++) { 01533 index[ind] = ind; 01534 } 01535 lastQ = MVT::CloneView( X, index ); 01536 01537 // Add these views to the Q and C arrays. 01538 Q.push_back( lastQ ); 01539 C.push_back( B ); 01540 qcs.push_back( MVT::GetNumberVecs( *lastQ ) ); 01541 } 01542 01543 // Get a view of the current vector in X to orthogonalize. 01544 indX[0] = j; 01545 Xj = MVT::CloneViewNonConst( X, indX ); 01546 if (this->_hasOp) { 01547 MXj = MVT::CloneViewNonConst( *MX, indX ); 01548 } 01549 else { 01550 MXj = Xj; 01551 } 01552 01553 // Compute the initial Op-norms 01554 MVT::MvDot( *Xj, *MXj, oldDot ); 01555 01556 Teuchos::Array<Teuchos::RCP<MV> > MQ(Q.size()); 01557 // Define the product Q^T * (Op*X) 01558 for (int i=0; i<Q.size(); i++) { 01559 01560 // Get a view of the current serial dense matrix 01561 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j ); 01562 01563 // Multiply Q' with MX 01564 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*Xj,MXj,tempC); 01565 // Multiply by Q and subtract the result in Xj 01566 MVT::MvTimesMatAddMv( -ONE, *Q[i], tempC, ONE, *Xj ); 01567 01568 // Update MXj, with the least number of applications of Op as possible 01569 if (this->_hasOp) { 01570 if (xc <= qcs[i]) { 01571 OPT::Apply( *(this->_Op), *Xj, *MXj); 01572 } 01573 else { 01574 // this will possibly be used again below; don't delete it 01575 MQ[i] = MVT::Clone( *Q[i], qcs[i] ); 01576 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] ); 01577 MVT::MvTimesMatAddMv( -ONE, *MQ[i], tempC, ONE, *MXj ); 01578 } 01579 } 01580 } 01581 01582 // Do any additional steps of classical Gram-Schmidt orthogonalization 01583 for (int num_ortho_steps=1; num_ortho_steps < max_ortho_steps_; ++num_ortho_steps) { 01584 01585 for (int i=0; i<Q.size(); i++) { 01586 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j ); 01587 Teuchos::SerialDenseMatrix<int,ScalarType> C2( qcs[i], 1 ); 01588 01589 // Apply another step of classical Gram-Schmidt 01590 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*Xj,MXj,C2); 01591 tempC += C2; 01592 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, *Xj ); 01593 01594 // Update MXj, with the least number of applications of Op as possible 01595 if (this->_hasOp) { 01596 if (MQ[i].get()) { 01597 // MQ was allocated and computed above; use it 01598 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MXj ); 01599 } 01600 else if (xc <= qcs[i]) { 01601 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01602 OPT::Apply( *(this->_Op), *Xj, *MXj); 01603 } 01604 } 01605 } // for (int i=0; i<Q.size(); i++) 01606 01607 } // for (int num_ortho_steps=1; num_ortho_steps < max_ortho_steps_; ++num_ortho_steps) 01608 01609 // Check for linear dependence. 01610 if (SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*sing_tol_)) { 01611 dep_flg = true; 01612 } 01613 01614 // Normalize the new vector if it's not dependent 01615 if (!dep_flg) { 01616 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01617 01618 MVT::MvScale( *Xj, ONE/diag ); 01619 if (this->_hasOp) { 01620 // Update MXj. 01621 MVT::MvScale( *MXj, ONE/diag ); 01622 } 01623 01624 // Enter value on diagonal of B. 01625 (*B)(j,j) = diag; 01626 } 01627 else { 01628 // Create a random vector and orthogonalize it against all previous columns of Q. 01629 Teuchos::RCP<MV> tempXj = MVT::Clone( X, 1 ); 01630 Teuchos::RCP<MV> tempMXj; 01631 MVT::MvRandom( *tempXj ); 01632 if (this->_hasOp) { 01633 tempMXj = MVT::Clone( X, 1 ); 01634 OPT::Apply( *(this->_Op), *tempXj, *tempMXj ); 01635 } 01636 else { 01637 tempMXj = tempXj; 01638 } 01639 MVT::MvDot( *tempXj, *tempMXj, oldDot ); 01640 // 01641 for (int num_orth=0; num_orth<max_ortho_steps_; num_orth++) { 01642 01643 for (int i=0; i<Q.size(); i++) { 01644 Teuchos::SerialDenseMatrix<int,ScalarType> product( qcs[i], 1 ); 01645 01646 // Apply another step of classical Gram-Schmidt 01647 MatOrthoManager<ScalarType,MV,OP>::innerProd(*Q[i],*tempXj,tempMXj,product); 01648 MVT::MvTimesMatAddMv( -ONE, *Q[i], product, ONE, *tempXj ); 01649 01650 // Update MXj, with the least number of applications of Op as possible 01651 if (this->_hasOp) { 01652 if (MQ[i].get()) { 01653 // MQ was allocated and computed above; use it 01654 MVT::MvTimesMatAddMv( -ONE, *MQ[i], product, ONE, *tempMXj ); 01655 } 01656 else if (xc <= qcs[i]) { 01657 // MQ was not allocated and computed above; it was cheaper to use X before and it still is 01658 OPT::Apply( *(this->_Op), *tempXj, *tempMXj); 01659 } 01660 } 01661 } // for (int i=0; i<nq; i++) 01662 01663 } 01664 01665 // Compute the Op-norms after the correction step. 01666 MVT::MvDot( *tempXj, *tempMXj, newDot ); 01667 01668 // Copy vector into current column of Xj 01669 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) { 01670 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0])); 01671 01672 // Enter value on diagonal of B. 01673 (*B)(j,j) = ZERO; 01674 01675 // Copy vector into current column of _basisvecs 01676 MVT::MvAddMv( ONE/diag, *tempXj, ZERO, *tempXj, *Xj ); 01677 if (this->_hasOp) { 01678 MVT::MvAddMv( ONE/diag, *tempMXj, ZERO, *tempMXj, *MXj ); 01679 } 01680 } 01681 else { 01682 return j; 01683 } 01684 } // if (!dep_flg) 01685 01686 // Remove the vectors from array 01687 if (j > 0) { 01688 Q.resize( nq ); 01689 C.resize( nq ); 01690 qcs.resize( nq ); 01691 } 01692 01693 } // for (int j=0; j<xc; j++) 01694 01695 return xc; 01696 } 01697 01698 } // namespace Belos 01699 01700 #endif // BELOS_ICGS_ORTHOMANAGER_HPP 01701
1.7.6.1