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