Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziTraceMinDavidson.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Anasazi: Block Eigensolvers Package
4 //
5 // Copyright 2004 NTESS and the Anasazi contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
14 #ifndef ANASAZI_TRACEMIN_DAVIDSON_HPP
15 #define ANASAZI_TRACEMIN_DAVIDSON_HPP
16 
17 #include "AnasaziConfigDefs.hpp"
18 #include "AnasaziEigensolver.hpp"
22 #include "AnasaziTraceMinBase.hpp"
23 
24 #include "Teuchos_ScalarTraits.hpp"
27 #include "Teuchos_TimeMonitor.hpp"
28 
29 
30 namespace Anasazi {
31 namespace Experimental {
32 
47  template <class ScalarType, class MV, class OP>
48  class TraceMinDavidson : public TraceMinBase<ScalarType,MV,OP> {
49  public:
50 
61  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
64  Teuchos::ParameterList &params
65  );
66 
67  private:
68  //
69  // Convenience typedefs
70  //
74  typedef typename SCT::magnitudeType MagnitudeType;
75 
76  // TraceMin specific methods
77  void addToBasis(const Teuchos::RCP<const MV> Delta);
78 
79  void harmonicAddToBasis(const Teuchos::RCP<const MV> Delta);
80  };
81 
84  //
85  // Implementations
86  //
89 
90 
91 
93  // Constructor
94  template <class ScalarType, class MV, class OP>
98  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
101  Teuchos::ParameterList &params
102  ) :
103  TraceMinBase<ScalarType,MV,OP>(problem,sorter,printer,tester,ortho,params)
104  {
105  }
106 
107 
109  // 1. Project Delta so that V' M Delta = 0 and Q' M Delta = 0
110  // 2. Normalize Delta so that Delta' M Delta = I
111  // 3. Add Delta to the end of V: V = [V Delta]
112  // 4. Update KV and MV
113  template <class ScalarType, class MV, class OP>
115  {
116  // TODO: We should also test the row length and map, etc...
117  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*Delta) != this->blockSize_, std::invalid_argument,
118  "Anasazi::TraceMinDavidson::addToBasis(): Delta does not have blockSize_ columns");
119 
120  int rank;
121  // Vector of indices
122  std::vector<int> curind(this->curDim_), newind(this->blockSize_);
123  // Pointer to the meaningful parts of V, KV, and MV
124  Teuchos::RCP<MV> lclV, lclKV, lclMV;
125  // Holds the vectors we project against
126  Teuchos::Array< Teuchos::RCP<const MV> > projVecs = this->auxVecs_;
127 
128  // Get the existing parts of the basis and add them to the list of things we project against
129  for(int i=0; i<this->curDim_; i++)
130  curind[i] = i;
131  lclV = MVT::CloneViewNonConst(*this->V_,curind);
132  projVecs.push_back(lclV);
133 
134  // Get the new part of the basis (where we're going to insert Delta)
135  for (int i=0; i<this->blockSize_; ++i)
136  newind[i] = this->curDim_ + i;
137  lclV = MVT::CloneViewNonConst(*this->V_,newind);
138 
139  // Insert Delta at the end of V
140  MVT::SetBlock(*Delta,newind,*this->V_);
141  this->curDim_ += this->blockSize_;
142 
143  // Project out the components of Delta in the direction of V
144  if(this->hasM_)
145  {
146  // It is more efficient to provide the orthomanager with MV
147  Teuchos::Array< Teuchos::RCP<const MV> > MprojVecs = this->MauxVecs_;
148  lclMV = MVT::CloneViewNonConst(*this->MV_,curind);
149  MprojVecs.push_back(lclMV);
150 
151  // Compute M*Delta
152  lclMV = MVT::CloneViewNonConst(*this->MV_,newind);
153  {
154  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
155  Teuchos::TimeMonitor lcltimer( *this->timerMOp_ );
156  #endif
157  this->count_ApplyM_+= this->blockSize_;
158  OPT::Apply(*this->MOp_,*lclV,*lclMV);
159  }
160 
161  {
162  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
163  Teuchos::TimeMonitor lcltimer( *this->timerOrtho_ );
164  #endif
165 
166  // Project and normalize Delta
167  rank = this->orthman_->projectAndNormalizeMat(*lclV,projVecs,
168  Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)),
169  Teuchos::null,lclMV,MprojVecs);
170  }
171 
172  MprojVecs.pop_back();
173  }
174  else
175  {
176  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
177  Teuchos::TimeMonitor lcltimer( *this->timerOrtho_ );
178  #endif
179 
180  // Project and normalize Delta
181  rank = this->orthman_->projectAndNormalizeMat(*lclV,projVecs);
182  }
183 
184  projVecs.pop_back();
185 
186  TEUCHOS_TEST_FOR_EXCEPTION(rank != this->blockSize_,TraceMinBaseOrthoFailure,
187  "Anasazi::TraceMinDavidson::addToBasis(): Couldn't generate basis of full rank.");
188 
189  // Update KV
190  if(this->Op_ != Teuchos::null)
191  {
192  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
193  Teuchos::TimeMonitor lcltimer( *this->timerOp_ );
194  #endif
195  this->count_ApplyOp_+= this->blockSize_;
196 
197  lclKV = MVT::CloneViewNonConst(*this->KV_,newind);
198  OPT::Apply(*this->Op_,*lclV,*lclKV);
199  }
200  }
201 
202 
203 
205  // 1. Project Delta so that V' M Delta = 0 and Q' M Delta = 0
206  // 2. Normalize Delta so that Delta' M Delta = I
207  // 3. Add Delta to the end of V: V = [V Delta]
208  // 4. Update KV and MV
209  template <class ScalarType, class MV, class OP>
210  void TraceMinDavidson<ScalarType,MV,OP>::harmonicAddToBasis(Teuchos::RCP<const MV> Delta)
211  {
212  // TODO: We should also test the row length and map, etc...
213  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*Delta) != this->blockSize_, std::invalid_argument,
214  "Anasazi::TraceMinDavidson::addToBasis(): Delta does not have blockSize_ columns");
215 
216  int rank;
217  // Vector of indices
218  std::vector<int> curind(this->curDim_), newind(this->blockSize_);
219  // Pointer to the meaningful parts of V, KV, and MV
220  Teuchos::RCP<MV> lclV, lclKV, lclMV, projVecs, KprojVecs;
221  // Array of things we project against
222 
223  // Get the existing parts of the basis and add them to the list of things we project against
224  for(int i=0; i<this->curDim_; i++)
225  curind[i] = i;
226  projVecs = MVT::CloneViewNonConst(*this->V_,curind);
227 
228  if(this->Op_ != Teuchos::null)
229  {
230  lclKV = MVT::CloneViewNonConst(*this->KV_,curind);
231  KprojVecs = lclKV;
232  }
233 
234  // Get the new part of the basis (where we're going to insert Delta)
235  for (int i=0; i<this->blockSize_; ++i)
236  newind[i] = this->curDim_ + i;
237  lclV = MVT::CloneViewNonConst(*this->V_,newind);
238 
239  // Insert Delta at the end of V
240  MVT::SetBlock(*Delta,newind,*this->V_);
241  this->curDim_ += this->blockSize_;
242 
243  // Project the auxVecs out of Delta
244  if(this->auxVecs_.size() > 0)
245  this->orthman_->projectMat(*lclV,this->auxVecs_);
246 
247  // Update KV
248  if(this->Op_ != Teuchos::null)
249  {
250  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
251  Teuchos::TimeMonitor lcltimer( *this->timerOp_ );
252  #endif
253  this->count_ApplyOp_+= this->blockSize_;
254 
255  lclKV = MVT::CloneViewNonConst(*this->KV_,newind);
256  OPT::Apply(*this->Op_,*lclV,*lclKV);
257  }
258 
259  // Project out the components of Delta in the direction of V
260 
261  // gamma = KauxVecs' lclKV
262  int nauxVecs = MVT::GetNumberVecs(*projVecs);
264 
265  this->orthman_->innerProdMat(*KprojVecs,*lclKV,*gamma);
266 
267  // lclKV = lclKV - KauxVecs gamma
268  MVT::MvTimesMatAddMv(-this->ONE,*KprojVecs,*gamma,this->ONE,*lclKV);
269 
270  // lclV = lclV - auxVecs gamma
271  MVT::MvTimesMatAddMv(-this->ONE,*projVecs,*gamma,this->ONE,*lclV);
272 
273  // Normalize lclKV
275  rank = this->orthman_->normalizeMat(*lclKV,gamma2);
276 
277  // lclV = lclV/gamma
279  SDsolver.setMatrix(gamma2);
280  SDsolver.invert();
281  RCP<MV> tempMV = MVT::CloneCopy(*lclV);
282  MVT::MvTimesMatAddMv(this->ONE,*tempMV,*gamma2,this->ZERO,*lclV);
283 
284  TEUCHOS_TEST_FOR_EXCEPTION(rank != this->blockSize_,TraceMinBaseOrthoFailure,
285  "Anasazi::TraceMinDavidson::addToBasis(): Couldn't generate basis of full rank.");
286 
287  // Update MV
288  if(this->hasM_)
289  {
290  #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
291  Teuchos::TimeMonitor lcltimer( *this->timerMOp_ );
292  #endif
293  this->count_ApplyM_+= this->blockSize_;
294 
295  lclMV = MVT::CloneViewNonConst(*this->MV_,newind);
296  OPT::Apply(*this->MOp_,*lclV,*lclMV);
297  }
298  }
299 
300 }} // End of namespace Anasazi
301 
302 #endif
303 
304 // End of file AnasaziTraceMinDavidson.hpp
Abstract base class for trace minimization eigensolvers.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
This class defines the interface required by an eigensolver and status test class to compute solution...
Declaration of basic traits for the multivector type.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Virtual base class which defines basic traits for the operator type.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
Anasazi&#39;s templated virtual class for providing routines for orthogonalization and orthonormalization...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Output managers remove the need for the eigensolver to know any information about the required output...
TraceMinDavidson(const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
TraceMinBase constructor with eigenproblem, solver utilities, and parameter list of solver options...
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
void push_back(const value_type &x)
This class implements a TraceMin-Davidson iteration for solving symmetric generalized eigenvalue prob...
This is an abstract base class for the trace minimization eigensolvers.
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Common interface of stopping criteria for Anasazi&#39;s solvers.
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)