Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
AnasaziTraceMinBase.hpp
Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 //
00004 //                 Anasazi: Block Eigensolvers Package
00005 //                 Copyright (2004) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 // TODO: Modify code so R is not computed unless needed
00030 
00035 #ifndef ANASAZI_TRACEMIN_BASE_HPP
00036 #define ANASAZI_TRACEMIN_BASE_HPP
00037 
00038 #include "AnasaziBasicSort.hpp"
00039 #include "AnasaziConfigDefs.hpp"
00040 #include "AnasaziEigensolver.hpp"
00041 #include "AnasaziMatOrthoManager.hpp"
00042 #include "AnasaziMultiVecTraits.hpp"
00043 #include "AnasaziOperatorTraits.hpp"
00044 #include "AnasaziSaddleOperator.hpp"
00045 #include "AnasaziSaddleContainer.hpp"
00046 #include "AnasaziSolverUtils.hpp"
00047 #include "AnasaziTraceMinRitzOp.hpp"
00048 #include "AnasaziTraceMinTypes.hpp"
00049 #include "AnasaziTypes.hpp"
00050 
00051 #include "Teuchos_ParameterList.hpp"
00052 #include "Teuchos_ScalarTraits.hpp"
00053 #include "Teuchos_SerialDenseMatrix.hpp"
00054 #include "Teuchos_TimeMonitor.hpp"
00055 
00056 using Teuchos::RCP;
00057 using Teuchos::rcp;
00058 
00059 namespace Anasazi {
00070 namespace Experimental {
00071 
00073 
00074 
00079   template <class ScalarType, class MV>
00080   struct TraceMinBaseState {
00082     int curDim;
00087     RCP<const MV> V;
00089     RCP<const MV> KV;
00091     RCP<const MV> MopV;
00093     RCP<const MV> X; 
00095     RCP<const MV> KX; 
00097     RCP<const MV> MX;
00099     RCP<const MV> R;
00101     RCP<const std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > T;
00107     RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > KK;
00109     RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > RV;
00111     bool isOrtho;
00113     int NEV;
00115     ScalarType largestSafeShift;
00117     RCP< const std::vector<ScalarType> > ritzShifts;
00118     TraceMinBaseState() : curDim(0), V(Teuchos::null), KV(Teuchos::null), MopV(Teuchos::null),
00119                           X(Teuchos::null), KX(Teuchos::null), MX(Teuchos::null), R(Teuchos::null), 
00120                           T(Teuchos::null), KK(Teuchos::null), RV(Teuchos::null), isOrtho(false), 
00121                           NEV(0), largestSafeShift(Teuchos::ScalarTraits<ScalarType>::zero()), 
00122                           ritzShifts(Teuchos::null) {}
00123   };
00124 
00126 
00128 
00129 
00142   class TraceMinBaseInitFailure : public AnasaziError {public:
00143     TraceMinBaseInitFailure(const std::string& what_arg) : AnasaziError(what_arg)
00144     {}};
00145 
00149   class TraceMinBaseOrthoFailure : public AnasaziError {public:
00150     TraceMinBaseOrthoFailure(const std::string& what_arg) : AnasaziError(what_arg)
00151     {}};
00152   
00154 
00167   template <class ScalarType, class MV, class OP>
00168   class TraceMinBase : public Eigensolver<ScalarType,MV,OP> { 
00169   public:
00171 
00172     
00199     TraceMinBase( const RCP<Eigenproblem<ScalarType,MV,OP> >    &problem, 
00200                    const RCP<SortManager<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > &sorter,
00201                    const RCP<OutputManager<ScalarType> >         &printer,
00202                    const RCP<StatusTest<ScalarType,MV,OP> >      &tester,
00203                    const RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00204                    Teuchos::ParameterList &params 
00205                  );
00206     
00208     virtual ~TraceMinBase();
00210 
00211 
00213 
00214 
00238     void iterate();
00239 
00262     void initialize(TraceMinBaseState<ScalarType,MV> newstate);
00263 
00267     void initialize();
00268 
00284     bool isInitialized() const;
00285 
00294     TraceMinBaseState<ScalarType,MV> getState() const;
00295     
00297 
00298 
00300 
00301 
00303     int getNumIters() const;
00304 
00306     void resetNumIters();
00307 
00315     RCP<const MV> getRitzVectors();
00316 
00322     std::vector<Value<ScalarType> > getRitzValues();
00323 
00324 
00333     std::vector<int> getRitzIndex();
00334 
00335 
00341     std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getResNorms();
00342 
00343 
00349     std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRes2Norms();
00350 
00351 
00359     std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRitzRes2Norms();
00360 
00366     int getCurSubspaceDim() const;
00367 
00369     int getMaxSubspaceDim() const;
00370 
00372 
00373 
00375 
00376 
00378     void setStatusTest(RCP<StatusTest<ScalarType,MV,OP> > test);
00379 
00381     RCP<StatusTest<ScalarType,MV,OP> > getStatusTest() const;
00382 
00384     const Eigenproblem<ScalarType,MV,OP>& getProblem() const;
00385 
00395     void setBlockSize(int blockSize);
00396 
00398     int getBlockSize() const;
00399 
00417     void setAuxVecs(const Teuchos::Array<RCP<const MV> > &auxvecs);
00418 
00420     Teuchos::Array<RCP<const MV> > getAuxVecs() const;
00421 
00423 
00425 
00426 
00433     void setSize(int blockSize, int numBlocks);
00434 
00436 
00437 
00439     void currentStatus(std::ostream &os);
00440 
00442 
00443   protected:
00444     //
00445     // Convenience typedefs
00446     //
00447     typedef SolverUtils<ScalarType,MV,OP>      Utils;
00448     typedef MultiVecTraits<ScalarType,MV>      MVT;
00449     typedef MultiVecTraitsExt<ScalarType,MV>   MVText;
00450     typedef OperatorTraits<ScalarType,MV,OP>   OPT;
00451     typedef Teuchos::ScalarTraits<ScalarType>  SCT;
00452     typedef typename SCT::magnitudeType        MagnitudeType;
00453     const MagnitudeType ONE;  
00454     const MagnitudeType ZERO; 
00455     const MagnitudeType NANVAL;
00456     //
00457     // Classes inputed through constructor that define the eigenproblem to be solved.
00458     //
00459     const RCP<Eigenproblem<ScalarType,MV,OP> >     problem_;
00460     const RCP<SortManager<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > sm_;
00461     const RCP<OutputManager<ScalarType> >          om_;
00462     RCP<StatusTest<ScalarType,MV,OP> >             tester_;
00463     const RCP<MatOrthoManager<ScalarType,MV,OP> >  orthman_;
00464     //
00465     // Information obtained from the eigenproblem
00466     //
00467     RCP<const OP> Op_;
00468     RCP<const OP> MOp_;
00469     RCP<const OP> Prec_;
00470     bool hasM_;
00471     //
00472     // Internal timers
00473     // TODO: Fix the timers
00474     //
00475     RCP<Teuchos::Time> timerOp_, timerMOp_, timerSaddle_, timerSortEval_, timerDS_,
00476                                 timerLocal_, timerCompRes_, timerOrtho_, timerInit_;
00477     //
00478     // Internal structs
00479     // TODO: Fix the checks
00480     //
00481     struct CheckList {
00482       bool checkV, checkX, checkMX, 
00483            checkKX, checkQ, checkKK;
00484       CheckList() : checkV(false),checkX(false),
00485                     checkMX(false),checkKX(false),
00486                     checkQ(false),checkKK(false) {};
00487     };
00488     //
00489     // Internal methods
00490     //
00491     std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
00492 
00493     //
00494     // Counters
00495     //
00496     int count_ApplyOp_, count_ApplyM_;
00497 
00498     //
00499     // Algorithmic parameters.
00500     //
00501     // blockSize_ is the solver block size; it controls the number of vectors added to the basis on each iteration.
00502     int blockSize_;
00503     // numBlocks_ is the size of the allocated space for the basis, in blocks.
00504     int numBlocks_; 
00505     
00506     // 
00507     // Current solver state
00508     //
00509     // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
00510     // is capable of running; _initialize is controlled  by the initialize() member method
00511     // For the implications of the state of initialized_, please see documentation for initialize()
00512     bool initialized_;
00513     //
00514     // curDim_ reflects how much of the current basis is valid 
00515     // NOTE: 0 <= curDim_ <= blockSize_*numBlocks_
00516     // this also tells us how many of the values in theta_ are valid Ritz values
00517     int curDim_;
00518     //
00519     // State Multivecs
00520     // V is the current basis
00521     // X is the current set of eigenvectors
00522     // R is the residual
00523     RCP<MV> X_, KX_, MX_, KV_, MV_, R_, V_;
00524     //
00525     // Projected matrices
00526     //
00527     RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > KK_, ritzVecs_;
00528     // 
00529     // auxiliary vectors
00530     Teuchos::Array<RCP<const MV> > auxVecs_, MauxVecs_;
00531     int numAuxVecs_;
00532     //
00533     // Number of iterations that have been performed.
00534     int iter_;
00535     // 
00536     // Current eigenvalues, residual norms
00537     std::vector<MagnitudeType> theta_, Rnorms_, R2norms_;
00538     // 
00539     // are the residual norms current with the residual?
00540     bool Rnorms_current_, R2norms_current_;
00541 
00542     // TraceMin specific variables
00543     RCP<TraceMinRitzOp<ScalarType,MV,OP> > ritzOp_;    // Operator which incorporates Ritz shifts
00544     enum SaddleSolType saddleSolType_;                          // Type of saddle point solver being used
00545     bool previouslyLeveled_;                                    // True if the trace already leveled off
00546     MagnitudeType previousTrace_;                               // Trace of X'KX at the previous iteration
00547     bool posSafeToShift_, negSafeToShift_;                      // Whether there are multiple clusters
00548     MagnitudeType largestSafeShift_;                            // The largest shift considered to be safe - generally the biggest converged eigenvalue
00549     int NEV_;                                                   // Number of eigenvalues we seek - used in computation of trace
00550     std::vector<ScalarType> ritzShifts_;                        // Current Ritz shifts
00551 
00552     // This is only used if we're using the Schur complement solver
00553     RCP<MV> Z_;
00554 
00555     // User provided TraceMin parameters
00556     enum WhenToShiftType whenToShift_;                          // What triggers a Ritz shift
00557     enum HowToShiftType howToShift_;                            // Strategy for choosing the Ritz shifts
00558     bool useMultipleShifts_;                                    // Whether to use one Ritz shift or many
00559     bool considerClusters_;                                     // Don't shift if there is only one cluster
00560     bool projectAllVecs_;                                       // Use all vectors in projected minres, or just 1
00561     bool projectLockedVecs_;                                    // Project locked vectors too in minres; does nothing if projectAllVecs = false
00562     bool computeAllRes_;                                        // Compute all residuals, or just blocksize ones - helps make Ritz shifts safer
00563     bool useRHSR_;                                              // Use -R as the RHS of projected minres rather than AX
00564     MagnitudeType traceThresh_;
00565 
00566     // Variables that are only used if we're shifting when the residual gets small
00567     // TODO: These are currently unused
00568     std::string shiftNorm_;                                     // Measure 2-norm or M-norm of residual
00569     MagnitudeType shiftThresh_;                                 // How small must the residual be?
00570     bool useRelShiftThresh_;                                    // Are we scaling the shift threshold by the eigenvalues?
00571 
00572     // TraceMin specific functions
00573     // Returns the trace of KK = X'KX
00574     ScalarType getTrace() const;
00575     // Returns true if the change in trace is very small (or was very small at one point)
00576     // TODO: Check whether I want to return true if the trace WAS small
00577     bool traceLeveled();
00578     // Get the residuals of each cluster of eigenvalues
00579     // TODO: Figure out whether I want to use these for all eigenvalues or just the first
00580     std::vector<ScalarType> getClusterResids();
00581     // Computes the Ritz shifts, which is not a trivial task
00582     // TODO: Make it safer for indefinite matrices maybe?
00583     void computeRitzShifts(const std::vector<ScalarType>& clusterResids);
00584     // Compute the tolerance for the inner solves
00585     // TODO: Come back to this and make sure it does what I want it to
00586     std::vector<ScalarType> computeTol();
00587     // Solves a saddle point problem
00588     void solveSaddlePointProblem(RCP<MV> Delta);
00589     // Solves a saddle point problem by using unpreconditioned projected minres
00590     void solveSaddleProj(RCP<MV> Delta) const;
00591     // Solves a saddle point problem by using preconditioned projected...Krylov...something
00592     // TODO: Fix this.  Replace minres with gmres or something.
00593     void solveSaddleProjPrec(RCP<MV> Delta) const;
00594     // Solves a saddle point problem by explicitly forming the inexact Schur complement
00595     void solveSaddleSchur (RCP<MV> Delta) const;
00596     // Solves a saddle point problem with a block diagonal preconditioner
00597     void solveSaddleBDPrec (RCP<MV> Delta) const; 
00598     // Computes KK = X'KX
00599     void computeKK();
00600     // Computes the eigenpairs of KK
00601     void computeRitzPairs();
00602     // Computes X = V evecs
00603     void computeX();
00604     // Updates KX and MX without a matvec by MAGIC (or by multiplying KV and MV by evecs)
00605     void updateKXMX();
00606     // Updates the residual
00607     void updateResidual();
00608     // Adds Delta to the basis
00609     // TraceMin and TraceMinDavidson each do this differently, which is why this is pure virtual
00610     virtual void addToBasis(const RCP<const MV> Delta) =0;
00611   };
00612 
00615   //
00616   // Implementations
00617   //
00620 
00621 
00623   // Constructor
00624   // TODO: Allow the users to supply their own linear solver
00625   // TODO: Add additional checking for logic errors (like trying to use gmres with multiple shifts)
00626   template <class ScalarType, class MV, class OP>
00627   TraceMinBase<ScalarType,MV,OP>::TraceMinBase(
00628         const RCP<Eigenproblem<ScalarType,MV,OP> >    &problem, 
00629         const RCP<SortManager<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > &sorter,
00630         const RCP<OutputManager<ScalarType> >         &printer,
00631         const RCP<StatusTest<ScalarType,MV,OP> >      &tester,
00632         const RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00633         Teuchos::ParameterList &params
00634         ) :
00635     ONE(Teuchos::ScalarTraits<MagnitudeType>::one()),
00636     ZERO(Teuchos::ScalarTraits<MagnitudeType>::zero()),
00637     NANVAL(Teuchos::ScalarTraits<MagnitudeType>::nan()),
00638     // problem, tools
00639     problem_(problem), 
00640     sm_(sorter),
00641     om_(printer),
00642     tester_(tester),
00643     orthman_(ortho),
00644     // timers, counters
00645 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
00646     timerOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Operation Op*x")),
00647     timerMOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Operation M*x")),
00648     timerSaddle_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Solving saddle point problem")),
00649     timerSortEval_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Sorting eigenvalues")),
00650     timerDS_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Direct solve")),
00651     timerLocal_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Local update")),
00652     timerCompRes_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Computing residuals")),
00653     timerOrtho_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Orthogonalization")),
00654     timerInit_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Initialization")),
00655 #endif
00656     count_ApplyOp_(0),
00657     count_ApplyM_(0),
00658     // internal data
00659     blockSize_(0),
00660     numBlocks_(0),
00661     initialized_(false),
00662     curDim_(0),
00663     auxVecs_( Teuchos::Array<RCP<const MV> >(0) ), 
00664     MauxVecs_( Teuchos::Array<RCP<const MV> >(0) ),
00665     numAuxVecs_(0),
00666     iter_(0),
00667     Rnorms_current_(false),
00668     R2norms_current_(false),
00669     // TraceMin variables
00670     previouslyLeveled_(false),
00671     previousTrace_(ZERO),
00672     posSafeToShift_(false),
00673     negSafeToShift_(false),
00674     largestSafeShift_(ZERO),
00675     Z_(Teuchos::null)
00676   {     
00677     TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,std::invalid_argument,
00678                        "Anasazi::TraceMinBase::constructor: user passed null problem pointer.");
00679     TEUCHOS_TEST_FOR_EXCEPTION(sm_ == Teuchos::null,std::invalid_argument,
00680                        "Anasazi::TraceMinBase::constructor: user passed null sort manager pointer.");
00681     TEUCHOS_TEST_FOR_EXCEPTION(om_ == Teuchos::null,std::invalid_argument,
00682                        "Anasazi::TraceMinBase::constructor: user passed null output manager pointer.");
00683     TEUCHOS_TEST_FOR_EXCEPTION(tester_ == Teuchos::null,std::invalid_argument,
00684                        "Anasazi::TraceMinBase::constructor: user passed null status test pointer.");
00685     TEUCHOS_TEST_FOR_EXCEPTION(orthman_ == Teuchos::null,std::invalid_argument,
00686                        "Anasazi::TraceMinBase::constructor: user passed null orthogonalization manager pointer.");
00687     TEUCHOS_TEST_FOR_EXCEPTION(problem_->isHermitian() == false, std::invalid_argument,
00688                        "Anasazi::TraceMinBase::constructor: problem is not hermitian.");
00689 
00690     // get the problem operators
00691     Op_   = problem_->getOperator();
00692     MOp_  = problem_->getM();
00693     Prec_ = problem_->getPrec();
00694     hasM_ = (MOp_ != Teuchos::null);
00695 
00696     // Set the saddle point solver parameters
00697     saddleSolType_ = params.get("Saddle Solver Type", PROJECTED_KRYLOV_SOLVER);
00698     TEUCHOS_TEST_FOR_EXCEPTION(saddleSolType_ != PROJECTED_KRYLOV_SOLVER && saddleSolType_ != SCHUR_COMPLEMENT_SOLVER && saddleSolType_ != BD_PREC_MINRES, std::invalid_argument,
00699            "Anasazi::TraceMin::constructor: Invalid value for \"Saddle Solver Type\"; valid options are PROJECTED_KRYLOV_SOLVER, SCHUR_COMPLEMENT_SOLVER, and BD_PREC_MINRES.");    
00700 
00701     // Set the Ritz shift parameters
00702     whenToShift_ = params.get("When To Shift", ALWAYS_SHIFT);
00703     TEUCHOS_TEST_FOR_EXCEPTION(whenToShift_ != NEVER_SHIFT && whenToShift_ != SHIFT_WHEN_TRACE_LEVELS && whenToShift_ != SHIFT_WHEN_RESID_SMALL && whenToShift_ != ALWAYS_SHIFT, std::invalid_argument,
00704            "Anasazi::TraceMin::constructor: Invalid value for \"When To Shift\"; valid options are \"NEVER_SHIFT\", \"SHIFT_WHEN_TRACE_LEVELS\", \"SHIFT_WHEN_RESID_SMALL\", and \"ALWAYS_SHIFT\".");  
00705 
00706     traceThresh_ = params.get("Trace Threshold", 2e-2);
00707     TEUCHOS_TEST_FOR_EXCEPTION(traceThresh_ < 0, std::invalid_argument,
00708            "Anasazi::TraceMin::constructor: Invalid value for \"Trace Threshold\"; Must be positive.");  
00709 
00710     howToShift_ = params.get("How To Choose Shift", ADJUSTED_RITZ_SHIFT);
00711     TEUCHOS_TEST_FOR_EXCEPTION(howToShift_ != LARGEST_CONVERGED_SHIFT && howToShift_ != ADJUSTED_RITZ_SHIFT && howToShift_ != RITZ_VALUES_SHIFT, std::invalid_argument,
00712            "Anasazi::TraceMin::constructor: Invalid value for \"How To Choose Shift\"; valid options are \"LARGEST_CONVERGED_SHIFT\", \"ADJUSTED_RITZ_SHIFT\", \"RITZ_VALUES_SHIFT\".");
00713 
00714     useMultipleShifts_ = params.get("Use Multiple Shifts", true);
00715 
00716     shiftThresh_ = params.get("Shift Threshold", 1e-2);
00717     useRelShiftThresh_ = params.get("Relative Shift Threshold", true);
00718     shiftNorm_ = params.get("Shift Norm", "2");
00719     TEUCHOS_TEST_FOR_EXCEPTION(shiftNorm_ != "2" && shiftNorm_ != "M", std::invalid_argument,
00720            "Anasazi::TraceMin::constructor: Invalid value for \"Shift Norm\"; valid options are \"2\", \"M\".");
00721 
00722     considerClusters_ = params.get("Consider Clusters", true);
00723 
00724     projectAllVecs_ = params.get("Project All Vectors", true);
00725     projectLockedVecs_ = params.get("Project Locked Vectors", true);
00726     useRHSR_ = params.get("Use Residual as RHS", true);
00727     computeAllRes_ = params.get("Compute All Residuals", true);
00728 
00729     // set the block size and allocate data
00730     int bs = params.get("Block Size", problem_->getNEV());
00731     int nb = params.get("Num Blocks", 1);
00732     setSize(bs,nb);
00733 
00734     NEV_ = problem_->getNEV();
00735 
00736     // Create the Ritz shift operator
00737     ritzOp_ = rcp( new TraceMinRitzOp<ScalarType,MV,OP>(Op_,MOp_,Prec_) );
00738   
00739   // Set the maximum number of inner iterations
00740   int innerMaxIts = params.get("Maximum Krylov Iterations", 1000);
00741   ritzOp_->setMaxIts(innerMaxIts);
00742   }
00743 
00744 
00746   // Destructor
00747   template <class ScalarType, class MV, class OP>
00748   TraceMinBase<ScalarType,MV,OP>::~TraceMinBase() {}
00749 
00750 
00752   // Set the block size
00753   // This simply calls setSize(), modifying the block size while retaining the number of blocks.
00754   template <class ScalarType, class MV, class OP>
00755   void TraceMinBase<ScalarType,MV,OP>::setBlockSize (int blockSize) 
00756   {
00757     TEUCHOS_TEST_FOR_EXCEPTION(blockSize < 1, std::invalid_argument, "Anasazi::TraceMinBase::setSize(blocksize,numblocks): blocksize must be strictly positive.");
00758     setSize(blockSize,numBlocks_);
00759   }
00760 
00761 
00763   // Return the current auxiliary vectors
00764   template <class ScalarType, class MV, class OP>
00765   Teuchos::Array<RCP<const MV> > TraceMinBase<ScalarType,MV,OP>::getAuxVecs() const {
00766     return auxVecs_;
00767   }
00768 
00769 
00771   // return the current block size
00772   template <class ScalarType, class MV, class OP>
00773   int TraceMinBase<ScalarType,MV,OP>::getBlockSize() const {
00774     return(blockSize_); 
00775   }
00776 
00777 
00779   // return eigenproblem
00780   template <class ScalarType, class MV, class OP>
00781   const Eigenproblem<ScalarType,MV,OP>& TraceMinBase<ScalarType,MV,OP>::getProblem() const { 
00782     return(*problem_); 
00783   }
00784 
00785 
00787   // return max subspace dim
00788   template <class ScalarType, class MV, class OP>
00789   int TraceMinBase<ScalarType,MV,OP>::getMaxSubspaceDim() const {
00790     return blockSize_*numBlocks_;
00791   }
00792 
00794   // return current subspace dim
00795   template <class ScalarType, class MV, class OP>
00796   int TraceMinBase<ScalarType,MV,OP>::getCurSubspaceDim() const {
00797     if (!initialized_) return 0;
00798     return curDim_;
00799   }
00800 
00801 
00803   // return ritz residual 2-norms
00804   template <class ScalarType, class MV, class OP>
00805   std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> 
00806   TraceMinBase<ScalarType,MV,OP>::getRitzRes2Norms() {
00807     return getRes2Norms();
00808   }
00809 
00810 
00812   // return ritz index
00813   template <class ScalarType, class MV, class OP>
00814   std::vector<int> TraceMinBase<ScalarType,MV,OP>::getRitzIndex() {
00815     std::vector<int> ret(curDim_,0);
00816     return ret;
00817   }
00818 
00819 
00821   // return ritz values
00822   template <class ScalarType, class MV, class OP>
00823   std::vector<Value<ScalarType> > TraceMinBase<ScalarType,MV,OP>::getRitzValues() { 
00824     std::vector<Value<ScalarType> > ret(curDim_);
00825     for (int i=0; i<curDim_; ++i) {
00826       ret[i].realpart = theta_[i];
00827       ret[i].imagpart = ZERO;
00828     }
00829     return ret;
00830   }
00831 
00832 
00834   // return pointer to ritz vectors
00835   template <class ScalarType, class MV, class OP>
00836   RCP<const MV> TraceMinBase<ScalarType,MV,OP>::getRitzVectors() {
00837     return X_;
00838   }
00839 
00840 
00842   // reset number of iterations
00843   template <class ScalarType, class MV, class OP>
00844   void TraceMinBase<ScalarType,MV,OP>::resetNumIters() { 
00845     iter_=0; 
00846   }
00847 
00848 
00850   // return number of iterations
00851   template <class ScalarType, class MV, class OP>
00852   int TraceMinBase<ScalarType,MV,OP>::getNumIters() const { 
00853     return(iter_); 
00854   }
00855 
00856 
00858   // return state pointers
00859   template <class ScalarType, class MV, class OP>
00860   TraceMinBaseState<ScalarType,MV> TraceMinBase<ScalarType,MV,OP>::getState() const {
00861     TraceMinBaseState<ScalarType,MV> state;
00862     state.curDim = curDim_;
00863     state.V = V_;
00864     state.X = X_;
00865     if (Op_ != Teuchos::null) {
00866       state.KV = KV_;
00867       state.KX = KX_;
00868     }
00869     else {
00870       state.KV = Teuchos::null;
00871       state.KX = Teuchos::null;
00872     }
00873     if (hasM_) {
00874       state.MopV = MV_;
00875       state.MX = MX_;
00876     }
00877     else {
00878       state.MopV = Teuchos::null;
00879       state.MX = Teuchos::null;
00880     }
00881     state.R = R_;
00882     state.KK = KK_;
00883     state.RV = ritzVecs_;
00884     if (curDim_ > 0) {
00885       state.T = rcp(new std::vector<MagnitudeType>(theta_.begin(),theta_.begin()+curDim_) );
00886     } else {
00887       state.T = rcp(new std::vector<MagnitudeType>(0));
00888     }
00889     state.ritzShifts = rcp(new std::vector<MagnitudeType>(ritzShifts_.begin(),ritzShifts_.begin()+blockSize_) );
00890     state.isOrtho = true;
00891     state.largestSafeShift = largestSafeShift_;
00892 
00893     return state;
00894   }
00895 
00896 
00898   // Perform TraceMinBase iterations until the StatusTest tells us to stop.
00899   template <class ScalarType, class MV, class OP>
00900   void TraceMinBase<ScalarType,MV,OP>::iterate ()
00901   {
00902     //
00903     // Initialize solver state
00904     if (initialized_ == false) {
00905       initialize();
00906     }
00907 
00908     // as a data member, this would be redundant and require synchronization with
00909     // blockSize_ and numBlocks_; we'll use a constant here.
00910     const int searchDim = blockSize_*numBlocks_;
00911 
00912     // Update obtained from solving saddle point problem
00913     RCP<MV> Delta = MVT::Clone(*X_,blockSize_);
00914 
00916     // iterate until the status test tells us to stop.
00917     // also break if our basis is full
00918     while (tester_->checkStatus(this) != Passed && (numBlocks_ == 1 || curDim_ < searchDim)) {
00919 
00920       // Print information on current iteration
00921       if (om_->isVerbosity(Debug)) {
00922         currentStatus( om_->stream(Debug) );
00923       }
00924       else if (om_->isVerbosity(IterationDetails)) {
00925         currentStatus( om_->stream(IterationDetails) );
00926       }
00927 
00928       ++iter_;
00929 
00930       // Solve the saddle point problem
00931       solveSaddlePointProblem(Delta);
00932 
00933       // Insert Delta at the end of V
00934       addToBasis(Delta);
00935 
00936       if (om_->isVerbosity( Debug ) ) {
00937         // Check almost everything here
00938         CheckList chk;
00939         chk.checkV = true;
00940         om_->print( Debug, accuracyCheck(chk, ": after addToBasis(Delta)") );
00941       }
00942 
00943       // Compute KK := V'KV
00944       computeKK();
00945 
00946       if (om_->isVerbosity( Debug ) ) {
00947         // Check almost everything here
00948         CheckList chk;
00949         chk.checkKK = true;
00950         om_->print( Debug, accuracyCheck(chk, ": after computeKK()") );
00951       }
00952 
00953       // Compute the Ritz pairs
00954       computeRitzPairs();
00955 
00956       // Compute X := V RitzVecs
00957       computeX();
00958 
00959       if (om_->isVerbosity( Debug ) ) {
00960         // Check almost everything here
00961         CheckList chk;
00962         chk.checkX = true;
00963         om_->print( Debug, accuracyCheck(chk, ": after computeX()") );
00964       }
00965 
00966       // Compute KX := KV RitzVecs and MX := MV RitzVecs (if necessary)
00967       updateKXMX();
00968 
00969       if (om_->isVerbosity( Debug ) ) {
00970         // Check almost everything here
00971         CheckList chk;
00972         chk.checkKX = true;
00973         chk.checkMX = true;
00974         om_->print( Debug, accuracyCheck(chk, ": after updateKXMX()") );
00975       }
00976 
00977       // Update the residual vectors
00978       updateResidual(); 
00979     } // end while (statusTest == false)
00980 
00981   } // end of iterate()
00982 
00983 
00985   // initialize the solver with default state
00986   template <class ScalarType, class MV, class OP>
00987   void TraceMinBase<ScalarType,MV,OP>::initialize()
00988   {
00989     TraceMinBaseState<ScalarType,MV> empty;
00990     initialize(empty);
00991   }
00992 
00993 
00995   /* Initialize the state of the solver
00996    * 
00997    * POST-CONDITIONS:
00998    *
00999    * V_ is orthonormal, orthogonal to auxVecs_, for first curDim_ vectors
01000    * theta_ contains Ritz w.r.t. V_(1:curDim_)
01001    * X is Ritz vectors w.r.t. V_(1:curDim_)
01002    * KX = Op*X
01003    * MX = M*X if hasM_
01004    * R = KX - MX*diag(theta_)
01005    *
01006    */
01007   template <class ScalarType, class MV, class OP>
01008   void TraceMinBase<ScalarType,MV,OP>::initialize(TraceMinBaseState<ScalarType,MV> newstate)
01009   {
01010     // TODO: Check for bad input
01011     // NOTE: memory has been allocated by setBlockSize(). Use setBlock below; do not Clone
01012     // NOTE: Overall time spent in this routine is counted to timerInit_; portions will also be counted towards other primitives
01013 
01014 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
01015     Teuchos::TimeMonitor inittimer( *timerInit_ );
01016 #endif
01017 
01018     std::vector<int> bsind(blockSize_);
01019     for (int i=0; i<blockSize_; ++i) bsind[i] = i;
01020 
01021     // in TraceMinBase, V is primary
01022     // the order of dependence follows like so.
01023     // --init->               V,KK
01024     //    --ritz analysis->   theta,X  
01025     //       --op apply->     KX,MX  
01026     //          --compute->   R
01027     // 
01028     // if the user specifies all data for a level, we will accept it.
01029     // otherwise, we will generate the whole level, and all subsequent levels.
01030     //
01031     // the data members are ordered based on dependence, and the levels are
01032     // partitioned according to the amount of work required to produce the
01033     // items in a level.
01034     //
01035     // inconsistent multivectors widths and lengths will not be tolerated, and
01036     // will be treated with exceptions.
01037     //
01038     // for multivector pointers in newstate which point directly (as opposed to indirectly, via a view) to 
01039     // multivectors in the solver, the copy will not be affected.
01040 
01041     // set up V and KK: get them from newstate if user specified them
01042     // otherwise, set them manually
01043     RCP<MV> lclV, lclKV, lclMV;
01044     RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclKK, lclRV;
01045 
01046     // If V was supplied by the user...
01047     if (newstate.V != Teuchos::null) {
01048       om_->stream(Debug) << "Copying V from the user\n";
01049 
01050       TEUCHOS_TEST_FOR_EXCEPTION( MVText::GetGlobalLength(*newstate.V) != MVText::GetGlobalLength(*V_), std::invalid_argument, 
01051              "Anasazi::TraceMinBase::initialize(newstate): Vector length of V not correct." );
01052       TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim < blockSize_, std::invalid_argument, 
01053              "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be at least blockSize().");
01054       TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*numBlocks_, std::invalid_argument, 
01055              "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be less than getMaxSubspaceDim().");
01056 
01057       TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < newstate.curDim, std::invalid_argument, 
01058              "Anasazi::TraceMinBase::initialize(newstate): Multivector for basis in new state must be as large as specified state rank.");
01059 
01060       curDim_ = newstate.curDim;
01061       // pick an integral amount
01062       curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
01063 
01064       TEUCHOS_TEST_FOR_EXCEPTION( curDim_ != newstate.curDim, std::invalid_argument, 
01065              "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be a multiple of getBlockSize().");
01066 
01067       // put data in V
01068       std::vector<int> nevind(curDim_);
01069       for (int i=0; i<curDim_; ++i) nevind[i] = i;
01070       if (newstate.V != V_) {
01071         if (curDim_ < MVT::GetNumberVecs(*newstate.V)) {
01072           newstate.V = MVT::CloneView(*newstate.V,nevind);
01073         }
01074         MVT::SetBlock(*newstate.V,nevind,*V_);
01075       }
01076       lclV = MVT::CloneViewNonConst(*V_,nevind);
01077     } // end if user specified V
01078     else 
01079     {
01080       // get vectors from problem or generate something, projectAndNormalize
01081       RCP<const MV> ivec = problem_->getInitVec();
01082       TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::invalid_argument,
01083              "Anasazi::TraceMinBase::initialize(newstate): Eigenproblem did not specify initial vectors to clone from.");
01084       // clear newstate so we won't use any data from it below
01085       newstate.X       = Teuchos::null;
01086       newstate.MX      = Teuchos::null;
01087       newstate.KX      = Teuchos::null;
01088       newstate.R       = Teuchos::null;
01089       newstate.T       = Teuchos::null;
01090       newstate.RV      = Teuchos::null;
01091       newstate.KK      = Teuchos::null;
01092       newstate.KV      = Teuchos::null;
01093       newstate.MopV    = Teuchos::null;
01094       newstate.isOrtho = false;
01095 
01096       // If the user did not specify a current dimension, use that of the initial vectors they provided
01097       if(newstate.curDim > 0)
01098         curDim_ = newstate.curDim;
01099       else
01100         curDim_ = MVT::GetNumberVecs(*ivec);
01101 
01102       // pick the largest multiple of blockSize_
01103       curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
01104       if (curDim_ > blockSize_*numBlocks_) {
01105         // user specified too many vectors... truncate
01106         // this produces a full subspace, but that is okay
01107         curDim_ = blockSize_*numBlocks_;
01108       }
01109       bool userand = false;
01110       if (curDim_ == 0) {
01111         // we need at least blockSize_ vectors
01112         // use a random multivec: ignore everything from InitVec
01113         userand = true;
01114         curDim_ = blockSize_;
01115       }
01116 
01117       std::vector<int> nevind(curDim_);
01118       for (int i=0; i<curDim_; ++i) nevind[i] = i;  
01119 
01120       // get a pointer into V
01121       // lclV has curDim vectors
01122       //
01123       // get pointer to first curDim vectors in V_
01124       lclV = MVT::CloneViewNonConst(*V_,nevind);
01125       if (userand) 
01126       {
01127         // generate random vector data
01128         MVT::MvRandom(*lclV);
01129       }
01130       else 
01131       {
01132         if(newstate.curDim > 0)
01133         {
01134           if(MVT::GetNumberVecs(*newstate.V) > curDim_) {
01135             RCP<const MV> helperMV = MVT::CloneView(*newstate.V,nevind);
01136             MVT::SetBlock(*helperMV,nevind,*lclV);
01137           }
01138           // assign ivec to first part of lclV
01139           MVT::SetBlock(*newstate.V,nevind,*lclV);
01140         }
01141         else
01142         {
01143           if(MVT::GetNumberVecs(*ivec) > curDim_) {
01144             ivec = MVT::CloneView(*ivec,nevind);
01145           }
01146           // assign ivec to first part of lclV
01147           MVT::SetBlock(*ivec,nevind,*lclV);
01148         }
01149       }
01150     } // end if user did not specify V
01151 
01152     // A vector of relevant indices
01153     std::vector<int> dimind(curDim_);
01154     for (int i=0; i<curDim_; ++i) dimind[i] = i;
01155 
01156     // Compute MV if necessary
01157     if(hasM_ && newstate.MopV == Teuchos::null)
01158     {
01159       om_->stream(Debug) << "Computing MV\n";
01160 
01161 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
01162       Teuchos::TimeMonitor lcltimer( *timerMOp_ );
01163 #endif
01164       count_ApplyM_+= curDim_;
01165 
01166       newstate.isOrtho = false;
01167       // Get a pointer to the relevant parts of MV
01168       lclMV = MVT::CloneViewNonConst(*MV_,dimind);
01169       OPT::Apply(*MOp_,*lclV,*lclMV);
01170     }
01171     // Copy MV if necessary
01172     else if(hasM_)
01173     {
01174       om_->stream(Debug) << "Copying MV\n";
01175 
01176       TEUCHOS_TEST_FOR_EXCEPTION( MVText::GetGlobalLength(*newstate.MopV) != MVText::GetGlobalLength(*MV_), std::invalid_argument, 
01177              "Anasazi::TraceMinBase::initialize(newstate): Vector length of MV not correct." );
01178 
01179       TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.MopV) < curDim_, std::invalid_argument, 
01180              "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in MV not correct.");
01181 
01182       if(newstate.MopV != MV_) {
01183         if(curDim_ < MVT::GetNumberVecs(*newstate.MopV)) {
01184           newstate.MopV = MVT::CloneView(*newstate.MopV,dimind);
01185         }
01186         MVT::SetBlock(*newstate.MopV,dimind,*MV_);
01187       }
01188       lclMV = MVT::CloneViewNonConst(*MV_,dimind);
01189     }
01190     // There is no M, so there is no MV
01191     else
01192     {
01193       om_->stream(Debug) << "There is no MV\n";
01194 
01195       lclMV = lclV;
01196       MV_ = V_;
01197     }
01198 
01199     // Project and normalize so that V'V = I and Q'V=0
01200     if(newstate.isOrtho == false)
01201     {
01202       om_->stream(Debug) << "Project and normalize\n";
01203 
01204 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
01205       Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
01206 #endif
01207 
01208       // These things are now invalid
01209       newstate.X       = Teuchos::null;
01210       newstate.MX      = Teuchos::null;
01211       newstate.KX      = Teuchos::null;
01212       newstate.R       = Teuchos::null;
01213       newstate.T       = Teuchos::null;
01214       newstate.RV      = Teuchos::null;
01215       newstate.KK      = Teuchos::null;
01216       newstate.KV      = Teuchos::null;
01217 
01218       int rank;
01219       if(auxVecs_.size() > 0)
01220       {
01221         rank = orthman_->projectAndNormalizeMat(*lclV, auxVecs_, 
01222                Teuchos::tuple(RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), 
01223                Teuchos::null, lclMV, MauxVecs_);
01224       }
01225       else
01226       {
01227         rank = orthman_->normalizeMat(*lclV,Teuchos::null,lclMV);
01228       }
01229 
01230       TEUCHOS_TEST_FOR_EXCEPTION(rank != curDim_,TraceMinBaseInitFailure,
01231              "Anasazi::TraceMinBase::initialize(): Couldn't generate initial basis of full rank.");
01232     }
01233 
01234     // Compute KV
01235     if(Op_ != Teuchos::null && newstate.KV == Teuchos::null)
01236     {
01237       om_->stream(Debug) << "Computing KV\n";
01238 
01239 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
01240       Teuchos::TimeMonitor lcltimer( *timerOp_ );
01241 #endif
01242       count_ApplyOp_+= curDim_;
01243 
01244       // These things are now invalid
01245       newstate.X       = Teuchos::null;
01246       newstate.MX      = Teuchos::null;
01247       newstate.KX      = Teuchos::null;
01248       newstate.R       = Teuchos::null;
01249       newstate.T       = Teuchos::null;
01250       newstate.RV      = Teuchos::null;
01251       newstate.KK      = Teuchos::null;
01252       newstate.KV      = Teuchos::null;
01253 
01254       lclKV = MVT::CloneViewNonConst(*KV_,dimind);
01255       OPT::Apply(*Op_,*lclV,*lclKV);
01256     }
01257     // Copy KV
01258     else if(Op_ != Teuchos::null)
01259     {
01260       om_->stream(Debug) << "Copying MV\n";
01261 
01262       TEUCHOS_TEST_FOR_EXCEPTION( MVText::GetGlobalLength(*newstate.KV) != MVText::GetGlobalLength(*KV_), std::invalid_argument, 
01263              "Anasazi::TraceMinBase::initialize(newstate): Vector length of MV not correct." );
01264 
01265       TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.KV) < curDim_, std::invalid_argument, 
01266              "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in KV not correct.");
01267 
01268       if (newstate.KV != KV_) {
01269         if (curDim_ < MVT::GetNumberVecs(*newstate.KV)) {
01270           newstate.KV = MVT::CloneView(*newstate.KV,dimind);
01271         }
01272         MVT::SetBlock(*newstate.KV,dimind,*KV_);
01273       }
01274       lclKV = MVT::CloneViewNonConst(*KV_,dimind);
01275     }
01276     else
01277     {
01278       om_->stream(Debug) << "There is no KV\n";
01279 
01280       lclKV = lclV;
01281       KV_ = V_;
01282     }      
01283 
01284     // Compute KK
01285     if(newstate.KK == Teuchos::null)
01286     {
01287       om_->stream(Debug) << "Computing KK\n";
01288 
01289       // These things are now invalid
01290       newstate.X       = Teuchos::null;
01291       newstate.MX      = Teuchos::null;
01292       newstate.KX      = Teuchos::null;
01293       newstate.R       = Teuchos::null;
01294       newstate.T       = Teuchos::null;
01295       newstate.RV      = Teuchos::null;
01296 
01297       // Note: there used to be a bug here.
01298       // We can't just call computeKK because it only computes the new parts of KK
01299 
01300       // Get a pointer to the part of KK we're interested in
01301       lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
01302 
01303       // KK := V'KV
01304       MVT::MvTransMv(ONE,*lclV,*lclKV,*lclKK);
01305     }
01306     // Copy KK
01307     else
01308     {
01309       om_->stream(Debug) << "Copying KK\n";
01310 
01311       // check size of KK
01312       TEUCHOS_TEST_FOR_EXCEPTION( newstate.KK->numRows() < curDim_ || newstate.KK->numCols() < curDim_, std::invalid_argument, 
01313              "Anasazi::TraceMinBase::initialize(newstate): Projected matrix in new state must be as large as specified state rank.");
01314 
01315       // put data into KK_
01316       lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
01317       if (newstate.KK != KK_) {
01318         if (newstate.KK->numRows() > curDim_ || newstate.KK->numCols() > curDim_) {
01319           newstate.KK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.KK,curDim_,curDim_) );
01320         }
01321         lclKK->assign(*newstate.KK);
01322       }
01323     }
01324 
01325     // Compute Ritz pairs
01326     if(newstate.T == Teuchos::null || newstate.RV == Teuchos::null)
01327     {
01328       om_->stream(Debug) << "Computing Ritz pairs\n";
01329 
01330       // These things are now invalid
01331       newstate.X       = Teuchos::null;
01332       newstate.MX      = Teuchos::null;
01333       newstate.KX      = Teuchos::null;
01334       newstate.R       = Teuchos::null;
01335       newstate.T       = Teuchos::null;
01336       newstate.RV      = Teuchos::null;
01337 
01338       computeRitzPairs();
01339     }
01340     // Copy Ritz pairs
01341     else
01342     {
01343       om_->stream(Debug) << "Copying Ritz pairs\n";
01344 
01345       TEUCHOS_TEST_FOR_EXCEPTION((signed int)(newstate.T->size()) != curDim_,
01346              std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of T must be consistent with dimension of V.");
01347 
01348       TEUCHOS_TEST_FOR_EXCEPTION( newstate.RV->numRows() < curDim_ || newstate.RV->numCols() < curDim_, std::invalid_argument, 
01349              "Anasazi::TraceMinBase::initialize(newstate): Ritz vectors in new state must be as large as specified state rank.");
01350 
01351       std::copy(newstate.T->begin(),newstate.T->end(),theta_.begin());
01352 
01353       lclRV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
01354       if (newstate.RV != ritzVecs_) {
01355         if (newstate.RV->numRows() > curDim_ || newstate.RV->numCols() > curDim_) {
01356           newstate.RV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.RV,curDim_,curDim_) );
01357         }
01358         lclRV->assign(*newstate.RV);
01359       }
01360     }
01361 
01362     // Compute X
01363     if(newstate.X == Teuchos::null)
01364     {
01365       om_->stream(Debug) << "Computing X\n";
01366 
01367       // These things are now invalid
01368       newstate.MX      = Teuchos::null;
01369       newstate.KX      = Teuchos::null;
01370       newstate.R       = Teuchos::null;
01371 
01372       computeX();
01373     }
01374     // Copy X
01375     else
01376     {
01377       om_->stream(Debug) << "Copying X\n";
01378 
01379       if(computeAllRes_ == false)
01380       {
01381         TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != blockSize_ || MVText::GetGlobalLength(*newstate.X) != MVText::GetGlobalLength(*X_),
01382                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with block size and length of V.");
01383 
01384         if(newstate.X != X_) {
01385           MVT::SetBlock(*newstate.X,bsind,*X_);
01386         }
01387       }
01388       else
01389       {
01390         TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != curDim_ || MVText::GetGlobalLength(*newstate.X) != MVText::GetGlobalLength(*X_),
01391                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with current dimension and length of V.");
01392 
01393         if(newstate.X != X_) {
01394           MVT::SetBlock(*newstate.X,dimind,*X_);
01395         }
01396       }
01397     }
01398 
01399     // Compute KX and MX if necessary
01400     // TODO: These technically should be separate; it won't matter much in terms of running time though
01401     if((Op_ != Teuchos::null && newstate.KX == Teuchos::null) || (hasM_ && newstate.MX == Teuchos::null))
01402     {
01403       om_->stream(Debug) << "Computing KX and MX\n";
01404 
01405       // These things are now invalid
01406       newstate.R = Teuchos::null;
01407 
01408       updateKXMX();
01409     }
01410     // Copy KX and MX if necessary
01411     else
01412     {
01413       om_->stream(Debug) << "Copying KX and MX\n";
01414       if(Op_ != Teuchos::null)
01415       {
01416         if(computeAllRes_ == false)
01417         {
01418           TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != blockSize_ || MVText::GetGlobalLength(*newstate.KX) != MVText::GetGlobalLength(*KX_),
01419                  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with block size and length of V.");
01420 
01421           if(newstate.KX != KX_) {
01422             MVT::SetBlock(*newstate.KX,bsind,*KX_);
01423           }
01424         }
01425         else
01426         {
01427           TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != curDim_ || MVText::GetGlobalLength(*newstate.KX) != MVText::GetGlobalLength(*KX_),
01428                  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with current dimension and length of V.");
01429 
01430           if (newstate.KX != KX_) {
01431             MVT::SetBlock(*newstate.KX,dimind,*KX_);
01432           }
01433         }
01434       }
01435 
01436       if(hasM_)
01437       {
01438         if(computeAllRes_ == false)
01439         {
01440           TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != blockSize_ || MVText::GetGlobalLength(*newstate.MX) != MVText::GetGlobalLength(*MX_),
01441                  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with block size and length of V.");
01442 
01443           if (newstate.MX != MX_) {
01444             MVT::SetBlock(*newstate.MX,bsind,*MX_);
01445           }
01446         }
01447         else
01448         {
01449           TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != curDim_ || MVText::GetGlobalLength(*newstate.MX) != MVText::GetGlobalLength(*MX_),
01450                  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with current dimension and length of V.");
01451 
01452           if (newstate.MX != MX_) {
01453             MVT::SetBlock(*newstate.MX,dimind,*MX_);
01454           }
01455         }
01456       }
01457     }
01458 
01459     // Compute R
01460     if(newstate.R == Teuchos::null)
01461     {
01462       om_->stream(Debug) << "Computing R\n";
01463 
01464       updateResidual();
01465     }
01466     // Copy R
01467     else
01468     {
01469       om_->stream(Debug) << "Copying R\n";
01470 
01471       if(computeAllRes_ == false)
01472       {
01473         TEUCHOS_TEST_FOR_EXCEPTION(MVText::GetGlobalLength(*newstate.R) != MVText::GetGlobalLength(*X_),
01474                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
01475         TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != blockSize_,
01476                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least block size vectors." );
01477         if (newstate.R != R_) {
01478           MVT::SetBlock(*newstate.R,bsind,*R_);
01479         }
01480       }
01481       else
01482       {
01483         TEUCHOS_TEST_FOR_EXCEPTION(MVText::GetGlobalLength(*newstate.R) != MVText::GetGlobalLength(*X_),
01484                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
01485         TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != curDim_,
01486                std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least curDim vectors." );
01487         if (newstate.R != R_) {
01488           MVT::SetBlock(*newstate.R,dimind,*R_);
01489         }
01490       }
01491     }
01492 
01493     // R has been updated; mark the norms as out-of-date
01494     Rnorms_current_ = false;
01495     R2norms_current_ = false;
01496 
01497     // Set the largest safe shift
01498     largestSafeShift_ = newstate.largestSafeShift;
01499 
01500     // Copy over the Ritz shifts
01501     if(newstate.ritzShifts != Teuchos::null)
01502     {
01503       om_->stream(Debug) << "Copying Ritz shifts\n";
01504       std::copy(newstate.ritzShifts->begin(),newstate.ritzShifts->end(),ritzShifts_.begin());
01505     }
01506     else
01507     {
01508       om_->stream(Debug) << "Setting Ritz shifts to 0\n";
01509       for(size_t i=0; i<ritzShifts_.size(); i++)
01510         ritzShifts_[i] = ZERO;
01511     }
01512 
01513     for(size_t i=0; i<ritzShifts_.size(); i++)
01514       om_->stream(Debug) << "Ritz shifts[" << i << "] = " << ritzShifts_[i] << std::endl;
01515 
01516     // finally, we are initialized
01517     initialized_ = true;
01518 
01519     if (om_->isVerbosity( Debug ) ) {
01520       // Check almost everything here
01521       CheckList chk;
01522       chk.checkV = true;
01523       chk.checkX = true;
01524       chk.checkKX = true;
01525       chk.checkMX = true;
01526       chk.checkQ = true;
01527       chk.checkKK = true;
01528       om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
01529     }
01530 
01531     // Print information on current status
01532     if (om_->isVerbosity(Debug)) {
01533       currentStatus( om_->stream(Debug) );
01534     }
01535     else if (om_->isVerbosity(IterationDetails)) {
01536       currentStatus( om_->stream(IterationDetails) );
01537     }
01538   }
01539 
01540 
01542   // Return initialized state
01543   template <class ScalarType, class MV, class OP>
01544   bool TraceMinBase<ScalarType,MV,OP>::isInitialized() const { return initialized_; }
01545 
01546 
01548   // Set the block size and make necessary adjustments.
01549   template <class ScalarType, class MV, class OP>
01550   void TraceMinBase<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks) 
01551   {
01552     // This routine only allocates space; it doesn't not perform any computation
01553     // any change in size will invalidate the state of the solver.
01554 
01555     TEUCHOS_TEST_FOR_EXCEPTION(blockSize < 1, std::invalid_argument, "Anasazi::TraceMinBase::setSize(blocksize,numblocks): blocksize must be strictly positive.");
01556 
01557     if (blockSize == blockSize_ && numBlocks == numBlocks_) {
01558       // do nothing
01559       return;
01560     }
01561 
01562     blockSize_ = blockSize;
01563     numBlocks_ = numBlocks;
01564 
01565     RCP<const MV> tmp;
01566     // grab some Multivector to Clone
01567     // in practice, getInitVec() should always provide this, but it is possible to use a 
01568     // Eigenproblem with nothing in getInitVec() by manually initializing with initialize(); 
01569     // in case of that strange scenario, we will try to Clone from X_ first, then resort to getInitVec()
01570     if (X_ != Teuchos::null) { // this is equivalent to blockSize_ > 0
01571       tmp = X_;
01572     }
01573     else {
01574       tmp = problem_->getInitVec();
01575       TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
01576              "Anasazi::TraceMinBase::setSize(): eigenproblem did not specify initial vectors to clone from.");
01577     }
01578 
01579     TEUCHOS_TEST_FOR_EXCEPTION(numAuxVecs_+blockSize*static_cast<ptrdiff_t>(numBlocks) > MVText::GetGlobalLength(*tmp),std::invalid_argument,
01580            "Anasazi::TraceMinBase::setSize(): max subspace dimension and auxilliary subspace too large.  Potentially impossible orthogonality constraints.");
01581 
01582     // New subspace dimension
01583     int newsd = blockSize_*numBlocks_;
01584 
01586     // blockSize dependent
01587     //
01588     ritzShifts_.resize(blockSize_,ZERO);
01589     if(computeAllRes_ == false)
01590     {
01591       // grow/allocate vectors
01592       Rnorms_.resize(blockSize_,NANVAL);
01593       R2norms_.resize(blockSize_,NANVAL);
01594       //
01595       // clone multivectors off of tmp
01596       //
01597       // free current allocation first, to make room for new allocation
01598       X_ = Teuchos::null;
01599       KX_ = Teuchos::null;
01600       MX_ = Teuchos::null;
01601       R_ = Teuchos::null;
01602       V_ = Teuchos::null;
01603       KV_ = Teuchos::null;
01604       MV_ = Teuchos::null;
01605 
01606       om_->print(Debug," >> Allocating X_\n");
01607       X_ = MVT::Clone(*tmp,blockSize_);
01608       if(Op_ != Teuchos::null) {
01609         om_->print(Debug," >> Allocating KX_\n");
01610         KX_ = MVT::Clone(*tmp,blockSize_);
01611       }
01612       else {
01613         KX_ = X_;
01614       }
01615       if (hasM_) {
01616         om_->print(Debug," >> Allocating MX_\n");
01617         MX_ = MVT::Clone(*tmp,blockSize_);
01618       }
01619       else {
01620         MX_ = X_;
01621       }
01622       om_->print(Debug," >> Allocating R_\n");
01623       R_ = MVT::Clone(*tmp,blockSize_);
01624     }
01625     else
01626     {
01627       // grow/allocate vectors
01628       Rnorms_.resize(newsd,NANVAL);
01629       R2norms_.resize(newsd,NANVAL);
01630       //
01631       // clone multivectors off of tmp
01632       //
01633       // free current allocation first, to make room for new allocation
01634       X_ = Teuchos::null;
01635       KX_ = Teuchos::null;
01636       MX_ = Teuchos::null;
01637       R_ = Teuchos::null;
01638       V_ = Teuchos::null;
01639       KV_ = Teuchos::null;
01640       MV_ = Teuchos::null;
01641 
01642       om_->print(Debug," >> Allocating X_\n");
01643       X_ = MVT::Clone(*tmp,newsd);
01644       if(Op_ != Teuchos::null) {
01645         om_->print(Debug," >> Allocating KX_\n");
01646         KX_ = MVT::Clone(*tmp,newsd);
01647       }
01648       else {
01649         KX_ = X_;
01650       }
01651       if (hasM_) {
01652         om_->print(Debug," >> Allocating MX_\n");
01653         MX_ = MVT::Clone(*tmp,newsd);
01654       }
01655       else {
01656         MX_ = X_;
01657       }
01658       om_->print(Debug," >> Allocating R_\n");
01659       R_ = MVT::Clone(*tmp,newsd);
01660     }
01661 
01663     // blockSize*numBlocks dependent
01664     //
01665     theta_.resize(newsd,NANVAL);
01666     om_->print(Debug," >> Allocating V_\n");
01667     V_ = MVT::Clone(*tmp,newsd);
01668     KK_ = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd,newsd) );
01669     ritzVecs_ = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd,newsd) );
01670 
01671     if(Op_ != Teuchos::null) {
01672       om_->print(Debug," >> Allocating KV_\n");
01673       KV_ = MVT::Clone(*tmp,newsd);
01674     }
01675     else {
01676       KV_ = V_;
01677     }
01678     if (hasM_) {
01679       om_->print(Debug," >> Allocating MV_\n");
01680       MV_ = MVT::Clone(*tmp,newsd);
01681     }
01682     else {
01683       MV_ = V_;
01684     }
01685 
01686     om_->print(Debug," >> done allocating.\n");
01687 
01688     initialized_ = false;
01689     curDim_ = 0;
01690   }
01691 
01692 
01694   // Set the auxiliary vectors
01695   template <class ScalarType, class MV, class OP>
01696   void TraceMinBase<ScalarType,MV,OP>::setAuxVecs(const Teuchos::Array<RCP<const MV> > &auxvecs) {
01697     typedef typename Teuchos::Array<RCP<const MV> >::iterator tarcpmv;
01698 
01699     // set new auxiliary vectors
01700     auxVecs_ = auxvecs;
01701 
01702     if(hasM_)
01703       MauxVecs_.resize(0);
01704     numAuxVecs_ = 0;
01705 
01706     for (tarcpmv i=auxVecs_.begin(); i != auxVecs_.end(); ++i) {
01707       numAuxVecs_ += MVT::GetNumberVecs(**i);
01708 
01709       if(hasM_)
01710       {
01711 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
01712       Teuchos::TimeMonitor lcltimer( *timerMOp_ );
01713 #endif
01714         count_ApplyM_+= MVT::GetNumberVecs(**i);
01715 
01716         RCP<MV> helperMV = MVT::Clone(**i,MVT::GetNumberVecs(**i));
01717         OPT::Apply(*MOp_,**i,*helperMV);
01718         MauxVecs_.push_back(helperMV);
01719       }
01720     }
01721 
01722     // If the solver has been initialized, V is not necessarily orthogonal to new auxiliary vectors
01723     if (numAuxVecs_ > 0 && initialized_) {
01724       initialized_ = false;
01725     }
01726 
01727     if (om_->isVerbosity( Debug ) ) {
01728       // Check almost everything here
01729       CheckList chk;
01730       chk.checkQ = true;
01731       om_->print( Debug, accuracyCheck(chk, ": after setAuxVecs()") );
01732     }
01733   }
01734 
01735 
01737   // compute/return residual M-norms
01738   template <class ScalarType, class MV, class OP>
01739   std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> 
01740   TraceMinBase<ScalarType,MV,OP>::getResNorms() {
01741     if (Rnorms_current_ == false) {
01742       // Update the residual norms
01743       if(computeAllRes_)
01744       {
01745         std::vector<int> curind(curDim_);
01746         for(int i=0; i<curDim_; i++)
01747           curind[i] = i;
01748 
01749         RCP<const MV> locR = MVT::CloneView(*R_,curind);
01750         std::vector<ScalarType> locNorms(curDim_);
01751         orthman_->norm(*locR,locNorms);
01752 
01753         for(int i=0; i<curDim_; i++)
01754           Rnorms_[i] = locNorms[i];
01755         for(int i=curDim_+1; i<blockSize_*numBlocks_; i++)
01756           Rnorms_[i] = NANVAL;
01757 
01758         Rnorms_current_ = true;
01759         locNorms.resize(blockSize_);
01760         return locNorms;
01761       }
01762       else
01763         orthman_->norm(*R_,Rnorms_);
01764       Rnorms_current_ = true;
01765     }
01766     else if(computeAllRes_)
01767     {
01768       std::vector<ScalarType> locNorms(blockSize_);
01769       for(int i=0; i<blockSize_; i++)
01770         locNorms[i] = Rnorms_[i];
01771       return locNorms;
01772     }
01773 
01774     return Rnorms_;
01775   }
01776 
01777   
01779   // compute/return residual 2-norms
01780   template <class ScalarType, class MV, class OP>
01781   std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> 
01782   TraceMinBase<ScalarType,MV,OP>::getRes2Norms() {
01783     if (R2norms_current_ == false) {
01784       // Update the residual 2-norms 
01785       if(computeAllRes_)
01786       {
01787         std::vector<int> curind(curDim_);
01788         for(int i=0; i<curDim_; i++)
01789           curind[i] = i;
01790 
01791         RCP<const MV> locR = MVT::CloneView(*R_,curind);
01792         std::vector<ScalarType> locNorms(curDim_);
01793         MVT::MvNorm(*locR,locNorms);
01794 
01795         for(int i=0; i<curDim_; i++)
01796         {
01797           R2norms_[i] = locNorms[i];
01798         }
01799         for(int i=curDim_+1; i<blockSize_*numBlocks_; i++)
01800           R2norms_[i] = NANVAL;
01801 
01802         R2norms_current_ = true;
01803         locNorms.resize(blockSize_);
01804         return locNorms;
01805       }
01806       else
01807         MVT::MvNorm(*R_,R2norms_);
01808       R2norms_current_ = true;
01809     }
01810     else if(computeAllRes_)
01811     {
01812       std::vector<ScalarType> locNorms(blockSize_);
01813       for(int i=0; i<blockSize_; i++)
01814         locNorms[i] = R2norms_[i];
01815       return locNorms;
01816     }
01817 
01818     return R2norms_;
01819   }
01820 
01822   // Set a new StatusTest for the solver.
01823   template <class ScalarType, class MV, class OP>
01824   void TraceMinBase<ScalarType,MV,OP>::setStatusTest(RCP<StatusTest<ScalarType,MV,OP> > test) {
01825     TEUCHOS_TEST_FOR_EXCEPTION(test == Teuchos::null,std::invalid_argument,
01826         "Anasazi::TraceMinBase::setStatusTest() was passed a null StatusTest.");
01827     tester_ = test;
01828   }
01829 
01831   // Get the current StatusTest used by the solver.
01832   template <class ScalarType, class MV, class OP>
01833   RCP<StatusTest<ScalarType,MV,OP> > TraceMinBase<ScalarType,MV,OP>::getStatusTest() const {
01834     return tester_;
01835   }
01836 
01838   // Print the current status of the solver
01839   template <class ScalarType, class MV, class OP>
01840   void 
01841   TraceMinBase<ScalarType,MV,OP>::currentStatus(std::ostream &os) 
01842   {
01843     using std::endl;
01844 
01845     os.setf(std::ios::scientific, std::ios::floatfield);
01846     os.precision(6);
01847     os <<endl;
01848     os <<"================================================================================" << endl;
01849     os << endl;
01850     os <<"                          TraceMinBase Solver Status" << endl;
01851     os << endl;
01852     os <<"The solver is "<<(initialized_ ? "initialized." : "not initialized.") << endl;
01853     os <<"The number of iterations performed is " <<iter_<<endl;
01854     os <<"The block size is         " << blockSize_<<endl;
01855     os <<"The number of blocks is   " << numBlocks_<<endl;
01856     os <<"The current basis size is " << curDim_<<endl;
01857     os <<"The number of auxiliary vectors is "<< numAuxVecs_ << endl;
01858     os <<"The number of operations Op*x   is "<<count_ApplyOp_<<endl;
01859     os <<"The number of operations M*x    is "<<count_ApplyM_<<endl;
01860 
01861     os.setf(std::ios_base::right, std::ios_base::adjustfield);
01862 
01863     if (initialized_) {
01864       os << endl;
01865       os <<"CURRENT EIGENVALUE ESTIMATES             "<<endl;
01866       os << std::setw(20) << "Eigenvalue" 
01867          << std::setw(20) << "Residual(M)"
01868          << std::setw(20) << "Residual(2)"
01869          << endl;
01870       os <<"--------------------------------------------------------------------------------"<<endl;
01871       for (int i=0; i<blockSize_; ++i) {
01872         os << std::setw(20) << theta_[i];
01873         if (Rnorms_current_) os << std::setw(20) << Rnorms_[i];
01874         else os << std::setw(20) << "not current";
01875         if (R2norms_current_) os << std::setw(20) << R2norms_[i];
01876         else os << std::setw(20) << "not current";
01877         os << endl;
01878       }
01879     }
01880     os <<"================================================================================" << endl;
01881     os << endl;
01882   }
01883 
01885   template <class ScalarType, class MV, class OP>
01886   ScalarType TraceMinBase<ScalarType,MV,OP>::getTrace() const
01887   {
01888     ScalarType currentTrace = ZERO;
01889 
01890     int neigs = std::min(NEV_,curDim_);
01891     for(int i=0; i < neigs; i++)
01892       currentTrace += theta_[i];
01893 
01894     return currentTrace;
01895   }
01896 
01898   template <class ScalarType, class MV, class OP>
01899   bool TraceMinBase<ScalarType,MV,OP>::traceLeveled()
01900   {
01901     ScalarType ratioOfChange = traceThresh_;
01902 
01903     if(previouslyLeveled_)
01904     {
01905       om_->stream(Debug) << "The trace already leveled, so we're not going to check it again\n";
01906       return true;
01907     }
01908 
01909     ScalarType currentTrace = getTrace();
01910 
01911     om_->stream(Debug) << "The current trace is " << currentTrace << std::endl;
01912 
01913     // Compute the ratio of the change
01914     // We seek the point where the trace has leveled off
01915     // It should be reasonably safe to shift at this point
01916     if(previousTrace_ != ZERO)
01917     {
01918       om_->stream(Debug) << "The previous trace was " << previousTrace_ << std::endl;
01919 
01920       ratioOfChange = std::abs(previousTrace_-currentTrace)/std::abs(previousTrace_);
01921       om_->stream(Debug) << "The ratio of change is " << ratioOfChange << std::endl;
01922     }
01923 
01924     previousTrace_ = currentTrace;
01925 
01926     if(ratioOfChange < traceThresh_)
01927     {
01928       previouslyLeveled_ = true;
01929       return true;
01930     }
01931     return false;
01932   }
01933 
01935   // Compute the residual of each CLUSTER of eigenvalues
01936   // This is important for selecting the Ritz shifts
01937   template <class ScalarType, class MV, class OP>
01938   std::vector<ScalarType> TraceMinBase<ScalarType,MV,OP>::getClusterResids()
01939   {
01940     int nvecs;
01941 
01942     if(computeAllRes_)
01943       nvecs = curDim_;
01944     else
01945       nvecs = blockSize_;
01946 
01947     std::vector<ScalarType> clusterResids(nvecs);
01948     std::vector<int> clusterIndices;
01949     for(int i=0; i < nvecs; i++)
01950     {
01951       // test for cluster
01952       if(clusterIndices.empty() || (theta_[i-1] + R2norms_[i-1] >= theta_[i] - R2norms_[i]))
01953       {
01954         // Add to cluster
01955         if(!clusterIndices.empty())  om_->stream(Debug) << theta_[i-1] << " is in a cluster with " << theta_[i] << " because " << theta_[i-1] + R2norms_[i-1] << " >= " << theta_[i] - R2norms_[i] << std::endl;
01956         clusterIndices.push_back(i);
01957       }
01958       // Cluster completed
01959       else
01960       {
01961         om_->stream(Debug) << theta_[i-1] << " is NOT in a cluster with " << theta_[i] << " because " << theta_[i-1] + R2norms_[i-1] << " < " << theta_[i] - R2norms_[i] << std::endl;
01962         ScalarType totalRes = ZERO;
01963         for(size_t j=0; j < clusterIndices.size(); j++)
01964           totalRes += (R2norms_[clusterIndices[j]]*R2norms_[clusterIndices[j]]);
01965 
01966         // If the smallest magnitude value of this sign is in a cluster with the 
01967         // largest magnitude cluster of this sign, it is not safe for the smallest
01968         // eigenvalue to use a shift
01969         if(theta_[clusterIndices[0]] < 0 && theta_[i] < 0)
01970           negSafeToShift_ = true;
01971         else if(theta_[clusterIndices[0]] > 0 && theta_[i] > 0)
01972           posSafeToShift_ = true;
01973 
01974         for(size_t j=0; j < clusterIndices.size(); j++)
01975           clusterResids[clusterIndices[j]] = sqrt(totalRes);
01976 
01977         clusterIndices.clear();
01978         clusterIndices.push_back(i);
01979       }
01980     }
01981 
01982     // Handle last cluster
01983     ScalarType totalRes = ZERO;
01984     for(size_t j=0; j < clusterIndices.size(); j++)
01985       totalRes += R2norms_[clusterIndices[j]];
01986     for(size_t j=0; j < clusterIndices.size(); j++)
01987       clusterResids[clusterIndices[j]] = totalRes;
01988 
01989     return clusterResids;
01990   }
01991 
01993   // Compute the Ritz shifts based on the Ritz values and residuals
01994   // A note on shifting: if the matrix is indefinite, you NEED to use a large block size
01995   // TODO: resids[i] on its own is unsafe for the generalized EVP
01996   //       See "A Parallel Implementation of the Trace Minimization Eigensolver"
01997   //       by Eloy Romero and Jose E. Roman
01998   template <class ScalarType, class MV, class OP>
01999   void TraceMinBase<ScalarType,MV,OP>::computeRitzShifts(const std::vector<ScalarType>& clusterResids)
02000   {
02001     std::vector<ScalarType> thetaMag(theta_);
02002     bool traceHasLeveled = traceLeveled();
02003 
02004     // Compute the magnitude of the eigenvalues
02005     for(int i=0; i<blockSize_; i++)
02006       thetaMag[i] = std::abs(thetaMag[i]);
02007 
02008     // Test whether it is safe to shift
02009     // TODO: Add residual shift option
02010     // Note: If you choose single shift with an indefinite matrix, you're gonna have a bad time...
02011     // Note: This is not safe for indefinite matrices, and I don't even know that it CAN be made safe.
02012     //       There just isn't any theoretical reason it should work.
02013     // TODO: It feels like this should be different for TraceMin-Base; we should be able to use all eigenvalues in the current subspace to determine whether we have a cluster
02014     if(whenToShift_ == ALWAYS_SHIFT || (whenToShift_ == SHIFT_WHEN_TRACE_LEVELS && traceHasLeveled))
02015     {
02016       // Set the shift to the largest safe shift
02017       if(howToShift_ == LARGEST_CONVERGED_SHIFT)
02018       {
02019         for(int i=0; i<blockSize_; i++)
02020           ritzShifts_[i] = largestSafeShift_;
02021       }
02022       // Set the shifts to the Ritz values
02023       else if(howToShift_ == RITZ_VALUES_SHIFT)
02024       {
02025         ritzShifts_[0] = theta_[0];
02026 
02027         // If we're using mulitple shifts, set them to EACH Ritz value.
02028         // Otherwise, only use the smallest
02029         if(useMultipleShifts_)
02030         {
02031           for(int i=1; i<blockSize_; i++)
02032             ritzShifts_[i] = theta_[i];
02033         }
02034         else
02035         {
02036           for(int i=1; i<blockSize_; i++)
02037             ritzShifts_[i] = ritzShifts_[0];
02038         }
02039       }
02040       // Use Dr. Sameh's original shifting strategy
02041       else if(howToShift_ == ADJUSTED_RITZ_SHIFT)
02042       {
02043         om_->stream(Debug) << "\nSeeking a shift for theta[0]=" << thetaMag[0] << std::endl;
02044 
02045         // This is my adjustment.  If all eigenvalues are in a single cluster, it's probably a bad idea to shift the smallest one.
02046         // If all of your eigenvalues are in one cluster, it's either way to early to shift or your subspace is too small
02047         if((theta_[0] > 0 && posSafeToShift_) || (theta_[0] < 0 && negSafeToShift_) || considerClusters_ == false)
02048         {
02049           // Initialize with a conservative shift, either the biggest safe shift or the eigenvalue adjusted by its cluster's residual
02050           ritzShifts_[0] = std::max(largestSafeShift_,thetaMag[0]-clusterResids[0]);
02051 
02052           om_->stream(Debug) << "Initializing with a conservative shift, either the most positive converged eigenvalue (" 
02053                              << largestSafeShift_ << ") or the eigenvalue adjusted by the residual (" << thetaMag[0] << "-" 
02054                              << clusterResids[0] << ").\n";
02055 
02056           // If this eigenvalue is NOT in a cluster, do an aggressive shift
02057           if(R2norms_[0] == clusterResids[0])
02058           {
02059             ritzShifts_[0] = thetaMag[0];
02060             om_->stream(Debug) << "Since this eigenvalue is NOT in a cluster, we can use the eigenvalue itself as a shift: ritzShifts[0]=" << ritzShifts_[0] << std::endl;
02061           }
02062           else
02063             om_->stream(Debug) << "This eigenvalue is in a cluster, so it would not be safe to use the eigenvalue itself as a shift\n";
02064         }
02065         else
02066         {
02067           if(largestSafeShift_ > std::abs(ritzShifts_[0]))
02068           {
02069             om_->stream(Debug) << "Initializing with a conservative shift...the most positive converged eigenvalue: " << largestSafeShift_ << std::endl;
02070             ritzShifts_[0] = largestSafeShift_;
02071           }
02072           else
02073             om_->stream(Debug) << "Using the previous value of ritzShifts[0]=" << ritzShifts_[0];
02074             
02075         }
02076 
02077         om_->stream(Debug) << "ritzShifts[0]=" << ritzShifts_[0] << std::endl;
02078 
02079         if(useMultipleShifts_)
02080         {
02082           // Compute shifts for other eigenvalues
02083           for(int i=1; i < blockSize_; i++)
02084           {
02085             om_->stream(Debug) << "\nSeeking a shift for theta[" << i << "]=" << thetaMag[i] << std::endl;
02086 
02087             // If the previous shift was aggressive and we are not in a cluster, do an aggressive shift
02088             if(ritzShifts_[i-1] == thetaMag[i-1] && i < blockSize_-1 && thetaMag[i] < thetaMag[i+1] - clusterResids[i+1])
02089             {
02090               ritzShifts_[i] = thetaMag[i];
02091               om_->stream(Debug) << "Using an aggressive shift: ritzShifts_[" << i << "]=" << ritzShifts_[i] << std::endl;
02092             }
02093             else
02094             {
02095               if(ritzShifts_[0] > std::abs(ritzShifts_[i]))
02096               {
02097                 om_->stream(Debug) << "It was unsafe to use the aggressive shift.  Choose the shift used by theta[0]=" 
02098                                    << thetaMag[0] << ": ritzShifts[0]=" << ritzShifts_[0] << std::endl;
02099 
02100                 // Choose a conservative shift, that of the smallest positive eigenvalue
02101                 ritzShifts_[i] = ritzShifts_[0];
02102               }
02103               else
02104                 om_->stream(Debug) << "It was unsafe to use the aggressive shift.  We will use the shift from the previous iteration: " << ritzShifts_[i] << std::endl;                
02105 
02106               om_->stream(Debug) << "Check whether any less conservative shifts would work (such as the biggest eigenvalue outside of the cluster, namely theta[ell] < " 
02107                                  << thetaMag[i] << "-" << clusterResids[i] << " (" << thetaMag[i] - clusterResids[i] << ")\n";
02108 
02109               // If possible, choose a less conservative shift, that of the biggest eigenvalue outside of the cluster
02110               for(int ell=0; ell < i; ell++)
02111               {
02112                 if(thetaMag[ell] < thetaMag[i] - clusterResids[i])
02113                 {
02114                   ritzShifts_[i] = thetaMag[ell];
02115                   om_->stream(Debug) << "ritzShifts_[" << i << "]=" << ritzShifts_[i] << " is valid\n";
02116                 }
02117                 else
02118                   break;
02119               }
02120             } // end else
02121 
02122             om_->stream(Debug) << "ritzShifts[" << i << "]=" << ritzShifts_[i] << std::endl;
02123           } // end for
02124         } // end if(useMultipleShifts_)
02125         else
02126         {
02127           for(int i=1; i<blockSize_; i++)
02128             ritzShifts_[i] = ritzShifts_[0];          
02129         }
02130       } // end if(howToShift_ == "Adjusted Ritz Values")
02131     } // end if(whenToShift_ == "Always" || (whenToShift_ == "After Trace Levels" && traceHasLeveled))
02132 
02133     // Set the correct sign
02134     for(int i=0; i<blockSize_; i++)
02135     {
02136       if(theta_[i] < 0)
02137         ritzShifts_[i] = -ritzShifts_[i];
02138     }
02139   }
02140 
02142   template <class ScalarType, class MV, class OP>
02143   std::vector<ScalarType> TraceMinBase<ScalarType,MV,OP>::computeTol()
02144   {
02145     ScalarType temp1, temp2;
02146     int nvecs = ritzShifts_.size();
02147     std::vector<ScalarType> tolerances;
02148   
02149     for(int i=0; i < nvecs; i++)
02150     {
02151       if(std::abs(theta_[0]) != std::abs(ritzShifts_[i]))
02152         temp1 = std::abs(theta_[i]-ritzShifts_[i])/std::abs(std::abs(theta_[0])-std::abs(ritzShifts_[i]));
02153       else
02154         temp1 = ONE;
02155 
02156       temp2 = pow(2.0,-iter_);
02157 
02158       // TODO: The min and max tolerances should not be hard coded
02159       //       Neither should the maximum number of iterations
02160       tolerances.push_back(std::max(std::min(temp1,temp2),1e-8));
02161     }
02162   
02163     if(nvecs > 1)
02164       tolerances[nvecs-1] = tolerances[nvecs-2];
02165 
02166     return tolerances;
02167   }
02168 
02170   template <class ScalarType, class MV, class OP>
02171   void TraceMinBase<ScalarType,MV,OP>::solveSaddlePointProblem(RCP<MV> Delta)
02172   {
02173 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02174   Teuchos::TimeMonitor lcltimer( *timerSaddle_ );
02175 #endif
02176 
02177     // This case can arise when looking for the largest eigenpairs
02178     if(Op_ == Teuchos::null)
02179     {
02180       // dense solver
02181       Teuchos::SerialDenseSolver<int,ScalarType> My_Solver;  
02182 
02183       // Schur complement
02184       RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclS = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(blockSize_,blockSize_) );
02185 
02186       if(computeAllRes_)
02187       {
02188         // Get the valid indices of X
02189         std::vector<int> curind(blockSize_);
02190         for(int i=0; i<blockSize_; i++)
02191           curind[i] = i;
02192 
02193         // Get a view of MX
02194         RCP<const MV> lclMX = MVT::CloneView(*MX_, curind);
02195 
02196         // form S = X' M^2 X
02197         MVT::MvTransMv(ONE,*lclMX,*lclMX,*lclS);
02198 
02199         // compute the inverse of S
02200         My_Solver.setMatrix(lclS);
02201         My_Solver.invert();
02202 
02203         // Delta = X - MX/S
02204         RCP<const MV> lclX = MVT::CloneView(*X_, curind);
02205         MVT::Assign(*lclX,*Delta);
02206         MVT::MvTimesMatAddMv( -ONE, *lclMX, *lclS, ONE, *Delta);
02207       }
02208       else
02209       {
02210         // form S = X' M^2 X
02211         MVT::MvTransMv(ONE,*MX_,*MX_,*lclS);
02212 
02213         // compute the inverse of S
02214         My_Solver.setMatrix(lclS);
02215         My_Solver.invert();
02216 
02217         // Delta = X - MX/S
02218         MVT::Assign(*X_,*Delta);
02219         MVT::MvTimesMatAddMv( -ONE, *MX_, *lclS, ONE, *Delta);
02220       }
02221     }
02222     else
02223     {
02224       std::vector<int> order(curDim_);
02225       std::vector<ScalarType> tempvec(blockSize_);
02226       RCP<BasicSort<MagnitudeType> > sorter = rcp( new BasicSort<MagnitudeType>("SR") );
02227 
02228       // Stores the residual of each CLUSTER of eigenvalues
02229       std::vector<ScalarType> clusterResids;
02230 
02231       // Sort the eigenvalues in ascending order for the Ritz shift selection
02232       sorter->sort(theta_, Teuchos::rcpFromRef(order), curDim_);   // don't catch exception
02233 
02234       // Apply the same ordering to the residual norms
02235       getRes2Norms();
02236       for (int i=0; i<blockSize_; i++)
02237         tempvec[i] = R2norms_[order[i]];
02238       R2norms_ = tempvec;
02239 
02240       // Compute the residual of each CLUSTER of eigenvalues
02241       // This is important for selecting the Ritz shifts
02242       clusterResids = getClusterResids();
02243 
02244       // Sort the eigenvalues based on what the user wanted
02245       sm_->sort(theta_, Teuchos::rcpFromRef(order), blockSize_);
02246 
02247       // Apply the same ordering to the residual norms and cluster residuals
02248       for (int i=0; i<blockSize_; i++)
02249         tempvec[i] = R2norms_[order[i]];
02250       R2norms_ = tempvec;
02251 
02252       for (int i=0; i<blockSize_; i++)
02253         tempvec[i] = clusterResids[order[i]];
02254       clusterResids = tempvec;      
02255 
02256       // Compute the Ritz shifts
02257       computeRitzShifts(clusterResids);
02258 
02259       // Compute the tolerances for the inner solves
02260       std::vector<ScalarType> tolerances = computeTol();
02261 
02262       for(int i=0; i<blockSize_; i++)
02263       {
02264         om_->stream(IterationDetails) << "Choosing Ritz shifts...theta[" << i << "]=" 
02265             << theta_[i] << ", resids[" << i << "]=" << R2norms_[i] << ", clusterResids[" << i << "]=" << clusterResids[i]
02266             << ", ritzShifts[" << i << "]=" << ritzShifts_[i] << ", and tol[" << i << "]=" << tolerances[i] << std::endl;
02267       }
02268 
02269       // Set the Ritz shifts for the solver
02270       ritzOp_->setRitzShifts(ritzShifts_);
02271 
02272       // Set the inner stopping tolerance
02273       // This uses the Ritz values to determine when to stop
02274       ritzOp_->setInnerTol(tolerances);
02275 
02276       // Solve the saddle point problem
02277       if(saddleSolType_ == PROJECTED_KRYLOV_SOLVER)
02278       {
02279         if(Prec_ != Teuchos::null)
02280           solveSaddleProjPrec(Delta);
02281         else
02282           solveSaddleProj(Delta);
02283       }
02284       else if(saddleSolType_ == SCHUR_COMPLEMENT_SOLVER)
02285       {
02286         if(Z_ == Teuchos::null || MVT::GetNumberVecs(*Z_) != blockSize_)
02287         {
02288           // We do NOT want Z to be 0, because that could result in stagnation
02289           Z_ = MVT::Clone(*X_,blockSize_);
02290           MVT::MvRandom(*Z_);
02291         }
02292         solveSaddleSchur(Delta);
02293       }
02294       else if(saddleSolType_ == BD_PREC_MINRES)
02295       {
02296         solveSaddleBDPrec(Delta);
02297 //        Delta->describe(*(Teuchos::VerboseObjectBase::getDefaultOStream()),Teuchos::VERB_EXTREME);
02298       }
02299       else
02300         std::cout << "Invalid saddle solver type\n";
02301     }
02302   }
02303 
02305   // Solve the saddle point problem using projected minres
02306   // TODO: We should be able to choose KX or -R as RHS.
02307   template <class ScalarType, class MV, class OP>
02308   void TraceMinBase<ScalarType,MV,OP>::solveSaddleProj (RCP<MV> Delta) const
02309   {
02310     RCP<TraceMinProjRitzOp<ScalarType,MV,OP> > projOp;
02311 
02312     if(computeAllRes_)
02313     {
02314       // Get the valid indices of X
02315       std::vector<int> curind(blockSize_);
02316       for(int i=0; i<blockSize_; i++)
02317         curind[i] = i;
02318 
02319       RCP<const MV> locMX = MVT::CloneView(*MX_, curind);
02320 
02321       // We should really project out the converged eigenvectors too
02322       if(projectAllVecs_)
02323       {
02324         if(projectLockedVecs_ && numAuxVecs_ > 0)
02325           projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX,orthman_,auxVecs_) );
02326         else
02327           projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX,orthman_) );
02328       }
02329       else
02330         projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX) );
02331 
02332       // Remember, Delta0 must equal 0
02333       // This ensures B-orthogonality between Delta and X
02334       MVT::MvInit(*Delta);
02335 
02336       if(useRHSR_)
02337       {
02338         RCP<const MV> locR = MVT::CloneView(*R_, curind);
02339         projOp->ApplyInverse(*locR, *Delta);
02340       }
02341       else
02342       {
02343         RCP<const MV> locKX = MVT::CloneView(*KX_, curind);
02344         projOp->ApplyInverse(*locKX, *Delta);
02345       }
02346     }
02347     else
02348     {
02349       // We should really project out the converged eigenvectors too
02350       if(projectAllVecs_)
02351       {
02352         if(projectLockedVecs_ && numAuxVecs_ > 0)
02353           projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_,orthman_,auxVecs_) );
02354         else
02355           projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_,orthman_) );
02356       }
02357       else
02358         projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_) );
02359 
02360       // Remember, Delta0 must equal 0
02361       // This ensures B-orthogonality between Delta and X
02362       MVT::MvInit(*Delta);
02363 
02364       if(useRHSR_) {
02365         projOp->ApplyInverse(*R_, *Delta);
02366       }
02367       else {
02368         projOp->ApplyInverse(*KX_, *Delta);
02369       }
02370     }
02371   }
02372 
02374   // TODO: Fix preconditioning
02375   template <class ScalarType, class MV, class OP>
02376   void TraceMinBase<ScalarType,MV,OP>::solveSaddleProjPrec (RCP<MV> Delta) const
02377   {
02378     // If we don't have Belos installed, we can't use TraceMinProjRitzOpWithPrec
02379     // Of course, this problem will be detected in the constructor and an exception will be thrown
02380     // This is only here to make sure the code compiles properly
02381 #ifdef HAVE_ANASAZI_BELOS
02382     RCP<TraceMinProjRitzOpWithPrec<ScalarType,MV,OP> > projOp;
02383 
02384     if(computeAllRes_)
02385     {
02386       // Get the valid indices of X
02387       std::vector<int> curind(blockSize_);
02388       for(int i=0; i<blockSize_; i++)
02389         curind[i] = i;
02390 
02391       RCP<const MV> locMX = MVT::CloneView(*MX_, curind);
02392 
02393       // We should really project out the converged eigenvectors too
02394       if(projectAllVecs_)
02395       {
02396         if(projectLockedVecs_ && numAuxVecs_ > 0)
02397           projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX,orthman_,auxVecs_) );
02398         else
02399           projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX,orthman_) );
02400       }
02401       else
02402         projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX) );
02403 
02404       // Remember, Delta0 must equal 0
02405       // This ensures B-orthogonality between Delta and X
02406       MVT::MvInit(*Delta);
02407 
02408       if(useRHSR_)
02409       {
02410         RCP<const MV> locR = MVT::CloneView(*R_, curind);
02411         projOp->ApplyInverse(*locR, *Delta);
02412         MVT::MvScale(*Delta, -ONE);
02413       }
02414       else
02415       {
02416         RCP<const MV> locKX = MVT::CloneView(*KX_, curind);
02417         projOp->ApplyInverse(*locKX, *Delta);
02418       }
02419     }
02420     else
02421     {
02422       // We should really project out the converged eigenvectors too
02423       if(projectAllVecs_)
02424       {
02425         if(projectLockedVecs_ && numAuxVecs_ > 0)
02426           projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_,orthman_,auxVecs_) );
02427         else
02428           projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_,orthman_) );
02429       }
02430       else
02431         projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_) );
02432 
02433       // Remember, Delta0 must equal 0
02434       // This ensures B-orthogonality between Delta and X
02435       MVT::MvInit(*Delta);
02436 
02437       if(useRHSR_)
02438       {
02439         projOp->ApplyInverse(*R_, *Delta);
02440         MVT::MvScale(*Delta,-ONE);
02441       }
02442       else
02443         projOp->ApplyInverse(*KX_, *Delta);
02444     }
02445 #endif
02446   }
02447 
02449   // TODO: We can hold the Schur complement constant in later iterations
02450   // TODO: Make sure we're using the preconditioner correctly
02451   template <class ScalarType, class MV, class OP>
02452   void TraceMinBase<ScalarType,MV,OP>::solveSaddleSchur (RCP<MV> Delta) const 
02453   {
02454     // dense solver
02455     Teuchos::SerialDenseSolver<int,ScalarType> My_Solver;  
02456 
02457     RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclL;
02458     RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclS = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(blockSize_,blockSize_) );
02459 
02460     if(computeAllRes_)
02461     {
02462       // Get the valid indices of X
02463       std::vector<int> curind(blockSize_);
02464       for(int i=0; i<blockSize_; i++)
02465         curind[i] = i;
02466 
02467       // Z = K \ MX
02468       MVT::SetBlock(*X_,curind,*Z_);
02469       RCP<const MV> lclMX = MVT::CloneView(*MX_, curind);
02470       ritzOp_->ApplyInverse(*lclMX,*Z_);
02471 
02472       // form S = X' M Z
02473       MVT::MvTransMv(ONE,*Z_,*lclMX,*lclS);
02474 
02475       // solve S L = I
02476       My_Solver.setMatrix(lclS);
02477       My_Solver.invert();
02478       lclL = lclS;
02479 
02480       // Delta = X - Z L
02481       RCP<const MV> lclX = MVT::CloneView(*X_, curind);
02482       MVT::Assign(*lclX,*Delta);
02483       MVT::MvTimesMatAddMv( -ONE, *Z_, *lclL, ONE, *Delta);
02484     }
02485     else
02486     {
02487       // Z = K \ MX
02488       ritzOp_->ApplyInverse(*MX_,*Z_);
02489 
02490       // form S = X' M Z
02491       MVT::MvTransMv(ONE,*Z_,*MX_,*lclS);
02492 
02493       // solve S L = I
02494       My_Solver.setMatrix(lclS);
02495       My_Solver.invert();
02496       lclL = lclS;
02497 
02498       // Delta = X - Z L
02499       MVT::Assign(*X_,*Delta);
02500       MVT::MvTimesMatAddMv( -ONE, *Z_, *lclL, ONE, *Delta);
02501     }
02502   }
02503 
02504 
02506   // TODO: We can hold the Schur complement constant in later iterations
02507   template <class ScalarType, class MV, class OP>
02508   void TraceMinBase<ScalarType,MV,OP>::solveSaddleBDPrec (RCP<MV> Delta) const 
02509   {
02510     RCP<MV> locKX, locMX;
02511     if(computeAllRes_)
02512     {
02513       std::vector<int> curind(blockSize_);
02514       for(int i=0; i<blockSize_; i++)
02515         curind[i] = i;
02516       locKX = MVT::CloneViewNonConst(*KX_, curind);
02517       locMX = MVT::CloneViewNonConst(*MX_, curind);
02518     }
02519     else
02520     {
02521       locKX = KX_;
02522       locMX = MX_;
02523     }
02524 
02525     // Create the operator [A BX; X'B 0]
02526     RCP< SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> > > sadOp = rcp(new SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> >(ritzOp_,locMX));
02527 
02528     // Create the RHS [AX; 0]
02529     RCP< SaddleContainer<ScalarType, MV> > sadRHS = rcp(new SaddleContainer<ScalarType,MV>(locKX));
02530 
02531     // Create the solution vector [Delta; L]
02532     MVT::MvInit(*Delta);
02533     RCP< SaddleContainer<ScalarType, MV> > sadSol = rcp(new SaddleContainer<ScalarType,MV>(Delta));
02534 
02535     // Create a minres solver
02536     RCP<PseudoBlockMinres<ScalarType,SaddleContainer<ScalarType, MV>,SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> > > > sadSolver;
02537     if(Prec_ != Teuchos::null)
02538     {
02539       RCP< SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> > > sadPrec = rcp(new SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> >(ritzOp_,locMX,BD_PREC));
02540       sadSolver = rcp(new PseudoBlockMinres<ScalarType,SaddleContainer<ScalarType, MV>,SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> > >(sadOp, sadPrec));
02541     }
02542     else {
02543       sadSolver = rcp(new PseudoBlockMinres<ScalarType,SaddleContainer<ScalarType, MV>,SaddleOperator<ScalarType, MV, TraceMinRitzOp<ScalarType,MV,OP> > >(sadOp));
02544     }
02545 
02546     // Set the tolerance for the minres solver
02547     std::vector<ScalarType> tol;
02548     ritzOp_->getInnerTol(tol);
02549     sadSolver->setTol(tol);
02550 
02551     // Set the maximum number of iterations
02552     sadSolver->setMaxIter(ritzOp_->getMaxIts());
02553 
02554     // Set the solution vector
02555     sadSolver->setSol(sadSol);
02556 
02557     // Set the RHS
02558     sadSolver->setRHS(sadRHS);
02559 
02560     // Solve the saddle point problem
02561     sadSolver->solve();
02562   }
02563 
02564 
02566   // Compute KK := V'KV
02567   // We only compute the NEW elements
02568   template <class ScalarType, class MV, class OP>
02569   void TraceMinBase<ScalarType,MV,OP>::computeKK()
02570   {
02571     // Get the valid indices of V
02572     std::vector<int> curind(curDim_);
02573     for(int i=0; i<curDim_; i++)
02574       curind[i] = i;
02575 
02576     // Get a pointer to the valid parts of V
02577     RCP<const MV> lclV = MVT::CloneView(*V_,curind);
02578 
02579     // Get the valid indices of KV
02580     curind.resize(blockSize_);
02581     for(int i=0; i<blockSize_; i++)
02582       curind[i] = curDim_-blockSize_+i;
02583     RCP<const MV> lclKV = MVT::CloneView(*KV_,curind);
02584 
02585     // Get a pointer to the valid part of KK
02586     RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > lclKK = 
02587         rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,blockSize_,0,curDim_-blockSize_) );
02588 
02589     // KK := V'KV
02590     MVT::MvTransMv(ONE,*lclV,*lclKV,*lclKK);
02591 
02592     // We only constructed the upper triangular part of the matrix, but that's okay because KK is symmetric!
02593     for(int r=0; r<curDim_; r++)
02594     {
02595       for(int c=0; c<r; c++)
02596       {
02597         (*KK_)(r,c) = (*KK_)(c,r);
02598       }
02599     }
02600   }
02601 
02603   // Compute the eigenpairs of KK, i.e. the Ritz pairs
02604   template <class ScalarType, class MV, class OP>
02605   void TraceMinBase<ScalarType,MV,OP>::computeRitzPairs()
02606   {
02607     // Get a pointer to the valid part of KK
02608     RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > lclKK = 
02609         rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
02610 
02611     // Get a pointer to the valid part of ritzVecs
02612     RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > lclRV = 
02613         rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
02614 
02615     // Compute Ritz pairs from KK
02616     {
02617 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02618       Teuchos::TimeMonitor lcltimer( *timerDS_ );
02619 #endif
02620       int rank = curDim_;
02621       Utils::directSolver(curDim_, *lclKK, Teuchos::null, *lclRV, theta_, rank, 10);
02622       // we want all ritz values back
02623       // TODO: This probably should not be an ortho failure
02624       TEUCHOS_TEST_FOR_EXCEPTION(rank != curDim_,TraceMinBaseOrthoFailure,
02625              "Anasazi::TraceMinBase::computeRitzPairs(): Failed to compute all eigenpairs of KK.");
02626     }
02627 
02628     // Sort ritz pairs
02629     {
02630 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02631       Teuchos::TimeMonitor lcltimer( *timerSortEval_ );
02632 #endif
02633 
02634       std::vector<int> order(curDim_);
02635       //
02636       // sort the first curDim_ values in theta_
02637       sm_->sort(theta_, Teuchos::rcpFromRef(order), curDim_);   // don't catch exception
02638       //
02639       // apply the same ordering to the primitive ritz vectors
02640       Utils::permuteVectors(order,*lclRV);
02641     }
02642   }
02643 
02645   // Compute X := V evecs
02646   template <class ScalarType, class MV, class OP>
02647   void TraceMinBase<ScalarType,MV,OP>::computeX()
02648   {
02649 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02650     Teuchos::TimeMonitor lcltimer( *timerLocal_ );
02651 #endif
02652 
02653     // Get the valid indices of V
02654     std::vector<int> curind(curDim_);
02655     for(int i=0; i<curDim_; i++)
02656       curind[i] = i;
02657 
02658     // Get a pointer to the valid parts of V
02659     RCP<const MV> lclV = MVT::CloneView(*V_,curind);
02660 
02661     if(computeAllRes_)
02662     {
02663       // Capture the relevant eigenvectors
02664       RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > relevantEvecs = 
02665           rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) ); 
02666 
02667       // X <- lclV*S
02668       RCP<MV> lclX = MVT::CloneViewNonConst(*X_,curind);
02669       MVT::MvTimesMatAddMv( ONE, *lclV, *relevantEvecs, ZERO, *lclX );
02670     }
02671     else
02672     {
02673       // Capture the relevant eigenvectors
02674       RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > relevantEvecs = 
02675           rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,blockSize_) ); 
02676 
02677       // X <- lclV*S
02678       MVT::MvTimesMatAddMv( ONE, *lclV, *relevantEvecs, ZERO, *X_ );
02679     }
02680   }
02681 
02683   // Compute KX := KV evecs and (if necessary) MX := MV evecs
02684   template <class ScalarType, class MV, class OP>
02685   void TraceMinBase<ScalarType,MV,OP>::updateKXMX()
02686   {
02687 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02688     Teuchos::TimeMonitor lcltimer( *timerLocal_ );
02689 #endif
02690 
02691     // Get the valid indices of V
02692     std::vector<int> curind(curDim_);
02693     for(int i=0; i<curDim_; i++)
02694       curind[i] = i;
02695 
02696     // Get pointers to the valid parts of V, KV, and MV (if necessary)
02697     RCP<const MV> lclV = MVT::CloneView(*V_,curind);
02698     RCP<const MV> lclKV = MVT::CloneView(*KV_,curind);
02699 
02700     if(computeAllRes_)
02701     {
02702       // Capture the relevant eigenvectors
02703       RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > relevantEvecs = 
02704           rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) ); 
02705 
02706       // Update KX and MX
02707       RCP<MV> lclKX = MVT::CloneViewNonConst(*KX_,curind);
02708       MVT::MvTimesMatAddMv( ONE, *lclKV, *relevantEvecs, ZERO, *lclKX );
02709       if(hasM_)
02710       {
02711         RCP<const MV> lclMV = MVT::CloneView(*MV_,curind);
02712         RCP<MV> lclMX = MVT::CloneViewNonConst(*MX_,curind);
02713         MVT::MvTimesMatAddMv( ONE, *lclMV, *relevantEvecs, ZERO, *lclMX );
02714       }
02715     }
02716     else
02717     {
02718       // Capture the relevant eigenvectors
02719       RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > relevantEvecs = 
02720           rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,blockSize_) ); 
02721 
02722       // Update KX and MX
02723       MVT::MvTimesMatAddMv( ONE, *lclKV, *relevantEvecs, ZERO, *KX_ );
02724       if(hasM_)
02725       {
02726         RCP<const MV> lclMV = MVT::CloneView(*MV_,curind);
02727         MVT::MvTimesMatAddMv( ONE, *lclMV, *relevantEvecs, ZERO, *MX_ );
02728       }
02729     }
02730   }
02731 
02733   // Update the residual R := KX - MX*T
02734   template <class ScalarType, class MV, class OP>
02735   void TraceMinBase<ScalarType,MV,OP>::updateResidual () {
02736 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
02737     Teuchos::TimeMonitor lcltimer( *timerCompRes_ );
02738 #endif
02739 
02740     if(computeAllRes_)
02741     {
02742       // Get the valid indices of X
02743       std::vector<int> curind(curDim_);
02744       for(int i=0; i<curDim_; i++)
02745         curind[i] = i;
02746 
02747       // Holds MX*T
02748       RCP<MV> MXT = MVT::CloneCopy(*MX_,curind);
02749 
02750       // Holds the relevant part of theta
02751       std::vector<ScalarType> locTheta(curDim_);
02752       for(int i=0; i<curDim_; i++)
02753         locTheta[i] = theta_[i];
02754 
02755       // Compute MX*T
02756       MVT::MvScale(*MXT,locTheta);
02757 
02758       // form R <- KX - MX*T
02759       RCP<const MV> locKX = MVT::CloneView(*KX_,curind);
02760       RCP<MV> locR = MVT::CloneViewNonConst(*R_,curind);
02761       MVT::MvAddMv(ONE,*locKX,-ONE,*MXT,*locR);
02762     }
02763     else
02764     {
02765       // Holds MX*T
02766       RCP<MV> MXT = MVT::CloneCopy(*MX_);
02767 
02768       // Holds the relevant part of theta
02769       std::vector<ScalarType> locTheta(blockSize_);
02770       for(int i=0; i<blockSize_; i++)
02771         locTheta[i] = theta_[i];
02772 
02773       // Compute MX*T
02774       MVT::MvScale(*MXT,locTheta);
02775 
02776       // form R <- KX - MX*T
02777       MVT::MvAddMv(ONE,*KX_,-ONE,*MXT,*R_);
02778     }
02779 
02780     // R has been updated; mark the norms as out-of-date
02781     Rnorms_current_ = false;
02782     R2norms_current_ = false;
02783   }
02784 
02785 
02787   // Check accuracy, orthogonality, and other debugging stuff
02788   // 
02789   // bools specify which tests we want to run (instead of running more than we actually care about)
02790   //
02791   // we don't bother checking the following because they are computed explicitly:
02792   //    H == Prec*R
02793   //   KH == K*H
02794   //
02795   // 
02796   // checkV : V orthonormal
02797   //          orthogonal to auxvecs
02798   // checkX : X orthonormal
02799   //          orthogonal to auxvecs
02800   // checkMX: check MX == M*X
02801   // checkKX: check KX == K*X
02802   // checkH : H orthonormal 
02803   //          orthogonal to V and H and auxvecs
02804   // checkMH: check MH == M*H
02805   // checkR : check R orthogonal to X
02806   // checkQ : check that auxiliary vectors are actually orthonormal
02807   // checkKK: check that KK is symmetric in memory 
02808   //
02809   // TODO: 
02810   //  add checkTheta 
02811   //
02812   template <class ScalarType, class MV, class OP>
02813   std::string TraceMinBase<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const 
02814   {
02815     using std::endl;
02816 
02817     std::stringstream os;
02818     os.precision(2);
02819     os.setf(std::ios::scientific, std::ios::floatfield);
02820 
02821     os << " Debugging checks: iteration " << iter_ << where << endl;
02822 
02823     // V and friends
02824     std::vector<int> lclind(curDim_);
02825     for (int i=0; i<curDim_; ++i) lclind[i] = i;
02826     RCP<const MV> lclV;
02827     if (initialized_) {
02828       lclV = MVT::CloneView(*V_,lclind);
02829     }
02830     if (chk.checkV && initialized_) {
02831       MagnitudeType err = orthman_->orthonormError(*lclV);
02832       os << " >> Error in V^H M V == I  : " << err << endl;
02833       for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
02834         err = orthman_->orthogError(*lclV,*auxVecs_[i]);
02835         os << " >> Error in V^H M Q[" << i << "] == 0 : " << err << endl;
02836       }
02837     }
02838 
02839     // X and friends
02840     RCP<const MV> lclX;
02841     if(initialized_)
02842     {
02843       if(computeAllRes_)
02844         lclX = MVT::CloneView(*X_,lclind);
02845       else
02846         lclX = X_;
02847     }
02848 
02849     if (chk.checkX && initialized_) {
02850       MagnitudeType err = orthman_->orthonormError(*lclX);
02851       os << " >> Error in X^H M X == I  : " << err << endl;
02852       for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
02853         err = orthman_->orthogError(*lclX,*auxVecs_[i]);
02854         os << " >> Error in X^H M Q[" << i << "] == 0 : " << err << endl;
02855       }
02856     }
02857     if (chk.checkMX && hasM_ && initialized_) {
02858       RCP<const MV> lclMX;
02859       if(computeAllRes_)
02860         lclMX = MVT::CloneView(*MX_,lclind);
02861       else
02862         lclMX = MX_;
02863 
02864       MagnitudeType err = Utils::errorEquality(*lclX, *lclMX, MOp_);
02865       os << " >> Error in MX == M*X     : " << err << endl;
02866     }
02867     if (Op_ != Teuchos::null && chk.checkKX && initialized_) {
02868       RCP<const MV> lclKX;
02869       if(computeAllRes_)
02870         lclKX = MVT::CloneView(*KX_,lclind);
02871       else
02872         lclKX = KX_;
02873 
02874       MagnitudeType err = Utils::errorEquality(*lclX, *lclKX, Op_);
02875       os << " >> Error in KX == K*X     : " << err << endl;
02876     }
02877 
02878     // KK
02879     if (chk.checkKK && initialized_) {
02880       Teuchos::SerialDenseMatrix<int,ScalarType> curKK(curDim_,curDim_);
02881       if(Op_ != Teuchos::null) {
02882         RCP<MV> lclKV = MVT::Clone(*V_,curDim_);
02883         OPT::Apply(*Op_,*lclV,*lclKV);
02884         MVT::MvTransMv(ONE,*lclV,*lclKV,curKK);
02885       }
02886       else {
02887         MVT::MvTransMv(ONE,*lclV,*lclV,curKK);
02888       }
02889       Teuchos::SerialDenseMatrix<int,ScalarType> subKK(Teuchos::View,*KK_,curDim_,curDim_);
02890       curKK -= subKK;
02891       os << " >> Error in V^H K V == KK : " << curKK.normFrobenius() << endl;
02892 
02893       Teuchos::SerialDenseMatrix<int,ScalarType> SDMerr(curDim_,curDim_);
02894       for (int j=0; j<curDim_; ++j) {
02895         for (int i=0; i<curDim_; ++i) {
02896           SDMerr(i,j) = subKK(i,j) - SCT::conjugate(subKK(j,i));
02897         }
02898       }
02899       os << " >> Error in KK - KK^H == 0 : " << SDMerr.normFrobenius() << endl;
02900     }
02901 
02902     // Q
02903     if (chk.checkQ) {
02904       for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
02905         MagnitudeType err = orthman_->orthonormError(*auxVecs_[i]);
02906         os << " >> Error in Q[" << i << "]^H M Q[" << i << "] == I : " << err << endl;
02907         for (Array_size_type j=i+1; j<auxVecs_.size(); ++j) {
02908           err = orthman_->orthogError(*auxVecs_[i],*auxVecs_[j]);
02909           os << " >> Error in Q[" << i << "]^H M Q[" << j << "] == 0 : " << err << endl;
02910         }
02911       }
02912     }
02913 
02914     os << endl;
02915 
02916     return os.str();
02917   }
02918   
02919 }} // End of namespace Anasazi
02920 
02921 #endif
02922 
02923 // End of file AnasaziTraceMinBase.hpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends