Belos Version of the Day
BelosSimpleOrthoManager.hpp
Go to the documentation of this file.
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 
00045 #ifndef __Belos_SimpleOrthoManager_hpp
00046 #define __Belos_SimpleOrthoManager_hpp
00047 
00048 #include <BelosConfigDefs.hpp>
00049 #include <BelosMultiVecTraits.hpp>
00050 #include <BelosOrthoManager.hpp>
00051 #include <BelosOutputManager.hpp>
00052 #include <Teuchos_ParameterList.hpp>
00053 #include <Teuchos_ParameterListAcceptorDefaultBase.hpp>
00054 #include <Teuchos_StandardCatchMacros.hpp>
00055 #include <Teuchos_TimeMonitor.hpp>
00056 
00057 namespace Belos {
00058 
00067   template<class Scalar, class MV>
00068   class SimpleOrthoManager : 
00069     public OrthoManager<Scalar, MV>, 
00070     public Teuchos::ParameterListAcceptorDefaultBase
00071   {
00072   public:
00073     typedef Scalar scalar_type;
00074     typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
00075     typedef Teuchos::SerialDenseMatrix<int, Scalar> mat_type;
00076     typedef Teuchos::RCP<Teuchos::SerialDenseMatrix<int, Scalar> > mat_ptr;
00077     
00078   private:
00079     typedef MultiVecTraits<Scalar, MV> MVT;
00080     typedef Teuchos::ScalarTraits<Scalar> STS;
00081     typedef Teuchos::ScalarTraits<magnitude_type> STM;
00082 
00084     std::string label_;
00086     Teuchos::RCP<OutputManager<Scalar> > outMan_;
00088     bool reorthogonalize_;
00090     bool useMgs_;
00092     mutable Teuchos::RCP<Teuchos::ParameterList> defaultParams_;
00093 
00094 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00095 
00096     Teuchos::RCP<Teuchos::Time> timerOrtho_;
00098     Teuchos::RCP<Teuchos::Time> timerProject_;
00100     Teuchos::RCP<Teuchos::Time> timerNormalize_;
00101 
00110     static Teuchos::RCP<Teuchos::Time>
00111     makeTimer (const std::string& prefix, 
00112          const std::string& timerName)
00113     {
00114       const std::string timerLabel = 
00115   prefix.empty() ? timerName : (prefix + ": " + timerName);
00116       return Teuchos::TimeMonitor::getNewTimer (timerLabel);
00117     }
00118 #endif // BELOS_TEUCHOS_TIME_MONITOR
00119     
00120   public:
00121 
00128     Teuchos::RCP<const Teuchos::ParameterList> 
00129     getValidParameters () const
00130     {
00131       using Teuchos::ParameterList;
00132       using Teuchos::parameterList;
00133       using Teuchos::RCP;
00134 
00135       const std::string defaultNormalizationMethod ("MGS");
00136       const bool defaultReorthogonalization = false;
00137 
00138       if (defaultParams_.is_null()) {
00139   RCP<ParameterList> params = parameterList ("Simple");
00140   params->set ("Normalization", defaultNormalizationMethod,
00141          "Which normalization method to use. Valid values are \"MGS\""
00142          " (for Modified Gram-Schmidt) and \"CGS\" (for Classical "
00143          "Gram-Schmidt).");
00144   params->set ("Reorthogonalization", defaultReorthogonalization,
00145          "Whether to perform one (unconditional) reorthogonalization "
00146          "pass.");
00147   defaultParams_ = params;
00148       }
00149       return defaultParams_;
00150     }
00151 
00161     static TEUCHOS_DEPRECATED Teuchos::RCP<const Teuchos::ParameterList> 
00162     getDefaultParameters ()
00163     {
00164       using Teuchos::ParameterList;
00165       using Teuchos::parameterList;
00166       using Teuchos::RCP;
00167 
00168       const std::string defaultNormalizationMethod ("MGS");
00169       const bool defaultReorthogonalization = false;
00170 
00171       RCP<ParameterList> params = parameterList ("Simple");
00172       params->set ("Normalization", defaultNormalizationMethod,
00173        "Which normalization method to use. Valid values are \"MGS\""
00174        " (for Modified Gram-Schmidt) and \"CGS\" (for Classical "
00175        "Gram-Schmidt).");
00176       params->set ("Reorthogonalization", defaultReorthogonalization,
00177        "Whether to perform one (unconditional) reorthogonalization "
00178        "pass.");
00179       return params;
00180     }
00181 
00187     Teuchos::RCP<const Teuchos::ParameterList> 
00188     getFastParameters ()
00189     {
00190       using Teuchos::ParameterList;
00191       using Teuchos::parameterList;
00192       using Teuchos::RCP;
00193       using Teuchos::rcp;
00194 
00195       const std::string fastNormalizationMethod ("CGS");
00196       const bool fastReorthogonalization = false;
00197 
00198       // Start with a clone of the default parameters.
00199       RCP<ParameterList> fastParams = parameterList (*getValidParameters());
00200       fastParams->set ("Normalization", fastNormalizationMethod);
00201       fastParams->set ("Reorthogonalization", fastReorthogonalization);
00202 
00203       return fastParams;
00204     }
00205 
00206     void 
00207     setParameterList (const Teuchos::RCP<Teuchos::ParameterList>& plist)
00208     {
00209       using Teuchos::ParameterList;
00210       using Teuchos::parameterList;
00211       using Teuchos::RCP;
00212       using Teuchos::Exceptions::InvalidParameter;
00213 
00214       RCP<const ParameterList> defaultParams = getValidParameters();
00215       RCP<ParameterList> params;
00216       if (plist.is_null()) {
00217   params = parameterList (*defaultParams);
00218       } else {
00219   params = plist;
00220   params->validateParametersAndSetDefaults (*defaultParams);
00221       }
00222       const std::string normalizeImpl = params->get<std::string>("Normalization");
00223       const bool reorthogonalize = params->get<bool>("Reorthogonalization");
00224 
00225       if (normalizeImpl == "MGS" || normalizeImpl == "Mgs" || normalizeImpl == "mgs") {
00226   useMgs_ = true;
00227   params->set ("Normalization", std::string ("MGS")); // Standardize.
00228       } else {
00229   useMgs_ = false;
00230   params->set ("Normalization", std::string ("CGS")); // Standardize.
00231       }
00232       reorthogonalize_ = reorthogonalize;
00233 
00234       setMyParamList (params);
00235     }
00236 
00247     SimpleOrthoManager (const Teuchos::RCP<OutputManager<Scalar> >& outMan,
00248       const std::string& label,
00249       const Teuchos::RCP<Teuchos::ParameterList>& params) :
00250       label_ (label),
00251       outMan_ (outMan)
00252     {
00253 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00254       timerOrtho_ = makeTimer (label, "All orthogonalization");
00255       timerProject_ = makeTimer (label, "Projection");
00256       timerNormalize_ = makeTimer (label, "Normalization");
00257 #endif // BELOS_TEUCHOS_TIME_MONITOR
00258 
00259       setParameterList (params);
00260       if (! outMan_.is_null()) {
00261   using std::endl;
00262   std::ostream& dbg = outMan_->stream(Debug);
00263   dbg << "Belos::SimpleOrthoManager constructor:" << endl
00264       << "-- Normalization method: " 
00265       << (useMgs_ ? "MGS" : "CGS") << endl
00266       << "-- Reorthogonalize (unconditionally)? " 
00267       << (reorthogonalize_ ? "Yes" : "No") << endl;
00268       }
00269     }
00270     
00272     virtual ~SimpleOrthoManager() {}
00273       
00274     void innerProd (const MV &X, const MV &Y, mat_type& Z ) const {
00275       MVT::MvTransMv (STS::one(), X, Y, Z);
00276     }
00277 
00278     void norm (const MV& X, std::vector<magnitude_type>& normVec) const {
00279       const int numCols = MVT::GetNumberVecs (X);
00280       // std::vector<T>::size_type is unsigned; int is signed.  Mixed
00281       // unsigned/signed comparisons trigger compiler warnings.
00282       if (normVec.size() < static_cast<size_t>(numCols))
00283   normVec.resize (numCols); // Resize normvec if necessary.
00284       MVT::MvNorm (X, normVec);
00285     }
00286 
00287     void 
00288     project (MV &X, 
00289        Teuchos::Array<mat_ptr> C,
00290        Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const
00291     {
00292 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00293       Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
00294       Teuchos::TimeMonitor timerMonitorProject(*timerProject_);
00295 #endif // BELOS_TEUCHOS_TIME_MONITOR
00296 
00297       allocateProjectionCoefficients (C, Q, X, true);
00298       rawProject (X, Q, C);
00299       if (reorthogonalize_) // Unconditional reorthogonalization
00300   {
00301     Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,Scalar> > > C2;
00302     allocateProjectionCoefficients (C2, Q, X, false);
00303     for (int k = 0; k < Q.size(); ++k)
00304       *C[k] += *C2[k];
00305   }
00306     }
00307 
00308     int 
00309     normalize (MV &X, mat_ptr B) const
00310     {
00311 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00312       Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
00313       Teuchos::TimeMonitor timerMonitorProject(*timerNormalize_);
00314 #endif // BELOS_TEUCHOS_TIME_MONITOR
00315 
00316       if (useMgs_)
00317   return normalizeMgs (X, B);
00318       else
00319   return normalizeCgs (X, B);
00320     }
00321 
00322     int 
00323     projectAndNormalize (MV &X, 
00324        Teuchos::Array<mat_ptr> C,
00325        mat_ptr B,
00326        Teuchos::ArrayView<Teuchos::RCP<const MV> > Q) const
00327     {
00328       // Don't need time monitors here: project() and normalize() have
00329       // their own.
00330       project (X, C, Q);
00331       return normalize (X, B);
00332     }
00333 
00334     magnitude_type 
00335     orthonormError(const MV &X) const 
00336     {
00337       const Scalar ONE = STS::one();
00338       const int ncols = MVT::GetNumberVecs(X);
00339       mat_type XTX (ncols, ncols);
00340       innerProd (X, X, XTX);
00341       for (int k = 0; k < ncols; ++k)
00342   XTX(k,k) -= ONE;
00343       return XTX.normFrobenius();
00344     }
00345 
00346     magnitude_type
00347     orthogError(const MV &X1, const MV &X2) const
00348     {
00349       const int ncols_X1 = MVT::GetNumberVecs (X1);
00350       const int ncols_X2 = MVT::GetNumberVecs (X2);
00351       mat_type X1_T_X2 (ncols_X1, ncols_X2);
00352       innerProd (X1, X2, X1_T_X2);
00353       return X1_T_X2.normFrobenius();
00354     }
00355 
00356     void setLabel (const std::string& label) { label_ = label; }
00357     const std::string& getLabel() const { return label_; }
00358 
00359   private:
00360 
00361     int
00362     normalizeMgs (MV &X, 
00363       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,Scalar> > B) const
00364     {
00365       using Teuchos::Range1D;
00366       using Teuchos::RCP;
00367       using Teuchos::rcp;
00368       using Teuchos::View;
00369 
00370       const int numCols = MVT::GetNumberVecs (X);
00371       if (numCols == 0)
00372   return 0;
00373 
00374       if (B.is_null())
00375   B = rcp (new mat_type (numCols, numCols));
00376       else if (B->numRows() != numCols || B->numCols() != numCols)
00377   B->shape (numCols, numCols);
00378 
00379       // Modified Gram-Schmidt orthogonalization
00380       std::vector<magnitude_type> normVec (1);
00381       for (int j = 0; j < numCols; ++j)
00382   {
00383     RCP<MV> X_cur = MVT::CloneViewNonConst (X, Range1D(j, j));
00384     MV& X_j = *X_cur;
00385     for (int i = 0; i < j; ++i)
00386       {
00387         RCP<const MV> X_prv = MVT::CloneView (X, Range1D(i, i));
00388         const MV& X_i = *X_prv;
00389         mat_type B_ij (View, *B, 1, 1, i, j);
00390         innerProd (X_i, X_j, B_ij);
00391         MVT::MvTimesMatAddMv (-STS::one(), X_i, B_ij, STS::one(), X_j);
00392         if (reorthogonalize_) // Unconditional reorthogonalization
00393     {       
00394       // innerProd() overwrites B(i,j), so save the
00395       // first-pass projection coefficient and update
00396       // B(i,j) afterwards.
00397       const Scalar B_ij_first = (*B)(i, j);
00398       innerProd (X_i, X_j, B_ij);
00399       MVT::MvTimesMatAddMv (-STS::one(), X_i, B_ij, STS::one(), X_j);
00400       (*B)(i, j) += B_ij_first;
00401     }
00402       }
00403     // Normalize column j of X
00404     norm (X_j, normVec);
00405     const magnitude_type theNorm = normVec[0];
00406     (*B)(j, j) = theNorm;
00407     if (normVec[0] != STM::zero())
00408       MVT::MvScale (X_j, STS::one() / theNorm);
00409     else
00410       return j; // break out early
00411   }
00412       return numCols; // full rank, as far as we know
00413     }
00414 
00415 
00416     int 
00417     normalizeCgs (MV &X, mat_ptr B) const
00418     {
00419       using Teuchos::Range1D;
00420       using Teuchos::RCP;
00421       using Teuchos::rcp;
00422       using Teuchos::View;
00423 
00424       const int numCols = MVT::GetNumberVecs (X);
00425       if (numCols == 0)
00426   return 0;
00427 
00428       if (B.is_null())
00429   B = rcp (new mat_type (numCols, numCols));
00430       else if (B->numRows() != numCols || B->numCols() != numCols)
00431   B->shape (numCols, numCols);
00432       mat_type& B_ref = *B;
00433 
00434       // Classical Gram-Schmidt orthogonalization 
00435       std::vector<magnitude_type> normVec (1);
00436 
00437       // Space for reorthogonalization
00438       mat_type B2 (numCols, numCols);
00439 
00440       // Do the first column first.
00441       {
00442   RCP<MV> X_cur = MVT::CloneViewNonConst (X, Range1D(0, 0));
00443   // Normalize column 0 of X
00444   norm (*X_cur, normVec);
00445   const magnitude_type theNorm = normVec[0];
00446   B_ref(0,0) = theNorm;
00447   if (theNorm != STM::zero())
00448     {
00449       const Scalar invNorm = STS::one() / theNorm;
00450       MVT::MvScale (*X_cur, invNorm);
00451     }
00452   else
00453     return 0; // break out early
00454       }
00455 
00456       // Orthogonalize the remaining columns of X
00457       for (int j = 1; j < numCols; ++j)
00458   {
00459     RCP<MV> X_cur = MVT::CloneViewNonConst (X, Range1D(j, j));
00460     RCP<const MV> X_prv = MVT::CloneView (X, Range1D(0, j-1));        
00461     mat_type B_prvcur (View, B_ref, j, 1, 0, j);
00462 
00463     // Project X_cur against X_prv (first pass)
00464     innerProd (*X_prv, *X_cur, B_prvcur);
00465     MVT::MvTimesMatAddMv (-STS::one(), *X_prv, B_prvcur, STS::one(), *X_cur);
00466     // Unconditional reorthogonalization: 
00467     // project X_cur against X_prv (second pass)
00468     if (reorthogonalize_) 
00469       {       
00470         mat_type B2_prvcur (View, B2, j, 1, 0, j);
00471         innerProd (*X_prv, *X_cur, B2_prvcur);
00472         MVT::MvTimesMatAddMv (-STS::one(), *X_prv, B2_prvcur, STS::one(), *X_cur);
00473         B_prvcur += B2_prvcur;
00474       }
00475     // Normalize column j of X
00476     norm (*X_cur, normVec);
00477     const magnitude_type theNorm = normVec[0];
00478     B_ref(j,j) = theNorm;
00479     if (theNorm != STM::zero())
00480       {
00481         const Scalar invNorm = STS::one() / theNorm;
00482         MVT::MvScale (*X_cur, invNorm);
00483       }
00484     else
00485       return j; // break out early
00486   }
00487       return numCols; // full rank, as far as we know
00488     }
00489 
00490 
00491     void
00492     allocateProjectionCoefficients (Teuchos::Array<mat_ptr>& C, 
00493             Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
00494             const MV& X,
00495             const bool attemptToRecycle = true) const
00496     {
00497       using Teuchos::rcp;
00498 
00499       const int num_Q_blocks = Q.size();
00500       const int ncols_X = MVT::GetNumberVecs (X);
00501       C.resize (num_Q_blocks);
00502       // # of block(s) that had to be (re)allocated (either allocated
00503       // freshly, or resized).
00504       int numAllocated = 0;
00505       if (attemptToRecycle)
00506   {
00507     for (int i = 0; i < num_Q_blocks; ++i) 
00508       {
00509         const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
00510         // Create a new C[i] if necessary, otherwise resize if
00511         // necessary, otherwise fill with zeros.
00512         if (C[i].is_null())
00513     {
00514       C[i] = rcp (new mat_type (ncols_Qi, ncols_X));
00515       numAllocated++;
00516     }
00517         else 
00518     {
00519       mat_type& Ci = *C[i];
00520       if (Ci.numRows() != ncols_Qi || Ci.numCols() != ncols_X)
00521         {
00522           Ci.shape (ncols_Qi, ncols_X);
00523           numAllocated++;
00524         }
00525       else
00526         Ci.putScalar (STS::zero());
00527     }
00528       }
00529   }
00530       else // Just allocate; don't try to check if we can recycle
00531   {
00532     for (int i = 0; i < num_Q_blocks; ++i) 
00533       {
00534         const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
00535         C[i] = rcp (new mat_type (ncols_Qi, ncols_X));
00536         numAllocated++;
00537       }
00538   }
00539       if (! outMan_.is_null())
00540   {
00541     using std::endl;
00542     std::ostream& dbg = outMan_->stream(Debug);
00543     dbg << "SimpleOrthoManager::allocateProjectionCoefficients: " 
00544         << "Allocated " << numAllocated << " blocks out of " 
00545         << num_Q_blocks << endl;
00546   }
00547     }
00548 
00549     void
00550     rawProject (MV& X, 
00551     Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
00552     Teuchos::ArrayView<mat_ptr> C) const
00553     { 
00554       // "Modified Gram-Schmidt" version of Block Gram-Schmidt.
00555       const int num_Q_blocks = Q.size();
00556       for (int i = 0; i < num_Q_blocks; ++i)
00557   {
00558     mat_type& Ci = *C[i];
00559     const MV& Qi = *Q[i];
00560     innerProd (Qi, X, Ci);
00561     MVT::MvTimesMatAddMv (-STS::one(), Qi, Ci, STS::one(), X);
00562   }
00563     }
00564 
00565   };
00566 } // namespace Belos
00567 
00568 #endif // __Belos_SimpleOrthoManager_hpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines