Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosCGIter.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_CG_ITER_HPP
11 #define BELOS_CG_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosCGIteration.hpp"
20 
21 #include "BelosLinearProblem.hpp"
22 #include "BelosOutputManager.hpp"
23 #include "BelosStatusTest.hpp"
25 #include "BelosOperatorTraits.hpp"
26 #include "BelosMultiVecTraits.hpp"
27 
30 #include "Teuchos_ScalarTraits.hpp"
32 #include "Teuchos_TimeMonitor.hpp"
33 
44 namespace Belos {
45 
46 template<class ScalarType, class MV, class OP>
47 class CGIter : virtual public CGIteration<ScalarType,MV,OP> {
48 
49  public:
50 
51  //
52  // Convenience typedefs
53  //
58 
60 
61 
68  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
71  Teuchos::ParameterList &params );
72 
74  virtual ~CGIter() {};
76 
77 
79 
80 
93  void iterate();
94 
110 
114  void initialize()
115  {
117  initializeCG(empty);
118  }
119 
128  state.R = R_;
129  state.P = P_;
130  state.AP = AP_;
131  state.Z = Z_;
132  return state;
133  }
134 
136 
137 
139 
140 
142  int getNumIters() const { return iter_; }
143 
145  void resetNumIters( int iter = 0 ) { iter_ = iter; }
146 
149  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> * norms ) const {
150  if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() == Belos::TwoNorm) {
151  (*norms)[0] = std::sqrt(Teuchos::ScalarTraits<ScalarType>::magnitude(rHr_));
152  return Teuchos::null;
153  } else
154  return R_;
155  }
156 
158 
160  Teuchos::RCP<MV> getCurrentUpdate() const { return Teuchos::null; }
161 
163 
165 
166 
168  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
169 
171  int getBlockSize() const { return 1; }
172 
174  void setBlockSize(int blockSize) {
175  TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
176  "Belos::CGIter::setBlockSize(): Cannot use a block size that is not one.");
177  }
178 
180  bool isInitialized() { return initialized_; }
181 
182 
184  void setDoCondEst(bool val) {
185  if (numEntriesForCondEst_) doCondEst_=val;
186  }
187 
190  // NOTE (mfh 30 Jul 2015) See note on getOffDiag() below.
191  // getDiag() didn't actually throw for me in that case, but why
192  // not be cautious?
193  typedef typename Teuchos::ArrayView<MagnitudeType>::size_type size_type;
194  if (static_cast<size_type> (iter_) >= diag_.size ()) {
195  return diag_ ();
196  } else {
197  return diag_ (0, iter_);
198  }
199  }
200 
203  // NOTE (mfh 30 Jul 2015) The implementation as I found it
204  // returned "offdiag(0,iter_)". This breaks (Teuchos throws in
205  // debug mode) when the maximum number of iterations has been
206  // reached, because iter_ == offdiag_.size() in that case. The
207  // new logic fixes this.
208  typedef typename Teuchos::ArrayView<MagnitudeType>::size_type size_type;
209  if (static_cast<size_type> (iter_) >= offdiag_.size ()) {
210  return offdiag_ ();
211  } else {
212  return offdiag_ (0, iter_);
213  }
214  }
215 
217 
218  private:
219 
220  //
221  // Internal methods
222  //
224  void setStateSize();
225 
226  //
227  // Classes inputed through constructor that define the linear problem to be solved.
228  //
233 
234  //
235  // Current solver state
236  //
237  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
238  // is capable of running; _initialize is controlled by the initialize() member method
239  // For the implications of the state of initialized_, please see documentation for initialize()
240  bool initialized_;
241 
242  // stateStorageInitialized_ specifies that the state storage has been initialized.
243  // This initialization may be postponed if the linear problem was generated without
244  // the right-hand side or solution vectors.
245  bool stateStorageInitialized_;
246 
247  // Current number of iterations performed.
248  int iter_;
249 
250  // Should the allreduce for convergence detection be merged with one of the inner products?
251  bool foldConvergenceDetectionIntoAllreduce_;
252 
253  // <r,r>
254  ScalarType rHr_;
255 
256  // Assert that the matrix is positive definite
257  bool assertPositiveDefiniteness_;
258 
259  // Tridiagonal system for condition estimation (if needed)
260  Teuchos::ArrayRCP<MagnitudeType> diag_, offdiag_;
261  int numEntriesForCondEst_;
262  bool doCondEst_;
263 
264 
265 
266  //
267  // State Storage
268  //
269  // Residual
270  Teuchos::RCP<MV> R_;
271  //
272  // Preconditioned residual
273  Teuchos::RCP<MV> Z_;
274  //
275  // Direction vector
276  Teuchos::RCP<MV> P_;
277  //
278  // Operator applied to direction vector
279  Teuchos::RCP<MV> AP_;
280 
281  Teuchos::RCP<MV> S_;
282 
283 };
284 
286  // Constructor.
287  template<class ScalarType, class MV, class OP>
289  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
292  Teuchos::ParameterList &params ):
293  lp_(problem),
294  om_(printer),
295  stest_(tester),
296  convTest_(convTester),
297  initialized_(false),
298  stateStorageInitialized_(false),
299  iter_(0),
300  assertPositiveDefiniteness_( params.get("Assert Positive Definiteness", true) ),
301  numEntriesForCondEst_(params.get("Max Size For Condest",0) ),
302  doCondEst_(false)
303  {
304  foldConvergenceDetectionIntoAllreduce_ = params.get<bool>("Fold Convergence Detection Into Allreduce",false);
305  }
306 
308  // Setup the state storage.
309  template <class ScalarType, class MV, class OP>
311  {
312  if (!stateStorageInitialized_) {
313 
314  // Check if there is any multivector to clone from.
315  Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
316  Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
317  if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
318  stateStorageInitialized_ = false;
319  return;
320  }
321  else {
322 
323  // Initialize the state storage
324  // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
325  if (R_ == Teuchos::null) {
326  // Get the multivector that is not null.
327  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
328  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
329  "Belos::CGIter::setStateSize(): linear problem does not specify multivectors to clone from.");
330  S_ = MVT::Clone( *tmp, 2 );
331  std::vector<int> index(1,0);
332  index[0] = 0;
333  R_ = MVT::CloneViewNonConst( *S_, index );
334  index[0] = 1;
335  Z_ = MVT::CloneViewNonConst( *S_, index );
336  P_ = MVT::Clone( *tmp, 1 );
337  AP_ = MVT::Clone( *tmp, 1 );
338 
339  }
340 
341  // Tracking information for condition number estimation
342  if(numEntriesForCondEst_ > 0) {
343  diag_.resize(numEntriesForCondEst_);
344  offdiag_.resize(numEntriesForCondEst_-1);
345  }
346 
347  // State storage has now been initialized.
348  stateStorageInitialized_ = true;
349  }
350  }
351  }
352 
353 
355  // Initialize this iteration object
356  template <class ScalarType, class MV, class OP>
358  {
359  // Initialize the state storage if it isn't already.
360  if (!stateStorageInitialized_)
361  setStateSize();
362 
363  TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
364  "Belos::CGIter::initialize(): Cannot initialize state storage!");
365 
366  // NOTE: In CGIter R_, the initial residual, is required!!!
367  //
368  std::string errstr("Belos::CGIter::initialize(): Specified multivectors must have a consistent length and width.");
369 
370  if (newstate.R != Teuchos::null) {
371 
372  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*R_),
373  std::invalid_argument, errstr );
374  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.R) != 1,
375  std::invalid_argument, errstr );
376 
377  // Copy basis vectors from newstate into V
378  if (newstate.R != R_) {
379  // copy over the initial residual (unpreconditioned).
380  MVT::Assign( *newstate.R, *R_ );
381  }
382 
383  // Compute initial direction vectors
384  // Initially, they are set to the preconditioned residuals
385  //
386  if ( lp_->getLeftPrec() != Teuchos::null ) {
387  lp_->applyLeftPrec( *R_, *Z_ );
388  if ( lp_->getRightPrec() != Teuchos::null ) {
389  Teuchos::RCP<MV> tmp = MVT::CloneCopy( *Z_ );
390  lp_->applyRightPrec( *tmp, *Z_ );
391  }
392  }
393  else if ( lp_->getRightPrec() != Teuchos::null ) {
394  lp_->applyRightPrec( *R_, *Z_ );
395  }
396  else {
397  MVT::Assign( *R_, *Z_ );
398  }
399  MVT::Assign( *Z_, *P_ );
400  }
401  else {
402 
403  TEUCHOS_TEST_FOR_EXCEPTION(newstate.R == Teuchos::null,std::invalid_argument,
404  "Belos::CGIter::initialize(): CGIterationState does not have initial residual.");
405  }
406 
407  // The solver is initialized
408  initialized_ = true;
409  }
410 
411 
413  // Iterate until the status test informs us we should stop.
414  template <class ScalarType, class MV, class OP>
416  {
417  //
418  // Allocate/initialize data structures
419  //
420  if (initialized_ == false) {
421  initialize();
422  }
423 
424  // Allocate memory for scalars.
425  std::vector<ScalarType> alpha(1);
426  std::vector<ScalarType> beta(1);
427  std::vector<ScalarType> rHz(1), rHz_old(1), pAp(1);
429 
430  // Create convenience variables for zero and one.
431  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
433 
434  // Scalars for condition estimation (if needed) - These will always use entry zero, for convenience
435  ScalarType pAp_old = one, beta_old = one, rHz_old2 = one;
436 
437  // Get the current solution vector.
438  Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
439 
440  // Check that the current solution vector only has one column.
441  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*cur_soln_vec) != 1, CGIterateFailure,
442  "Belos::CGIter::iterate(): current linear system has more than one vector!" );
443 
444  // Compute first <r,z> a.k.a. rHz
445  if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() == Belos::TwoNorm) {
446  MVT::MvTransMv( one, *R_, *S_, rHs );
447  rHr_ = rHs(0,0);
448  rHz[0] = rHs(0,1);
449  } else
450  MVT::MvDot( *R_, *Z_, rHz );
451 
453  // Iterate until the status test tells us to stop.
454  //
455  while (stest_->checkStatus(this) != Passed) {
456 
457  // Increment the iteration
458  iter_++;
459 
460  // Multiply the current direction vector by A and store in AP_
461  lp_->applyOp( *P_, *AP_ );
462 
463  // Compute alpha := <R_,Z_> / <P_,AP_>
464  MVT::MvDot( *P_, *AP_, pAp );
465  alpha[0] = rHz[0] / pAp[0];
466 
467  // Check that alpha is a positive number!
468  if(assertPositiveDefiniteness_) {
469  TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(alpha[0]) <= zero, CGPositiveDefiniteFailure,
470  "Belos::CGIter::iterate(): non-positive value for p^H*A*p encountered!" );
471  }
472 
473  //
474  // Update the solution vector x := x + alpha * P_
475  //
476  MVT::MvAddMv( one, *cur_soln_vec, alpha[0], *P_, *cur_soln_vec );
477  lp_->updateSolution();
478  //
479  // Save the denominator of beta before residual is updated [ old <R_, Z_> ]
480  //
481  rHz_old[0] = rHz[0];
482  //
483  // Compute the new residual R_ := R_ - alpha * AP_
484  //
485  MVT::MvAddMv( one, *R_, -alpha[0], *AP_, *R_ );
486  //
487  // Compute beta := [ new <R_, Z_> ] / [ old <R_, Z_> ],
488  // and the new direction vector p.
489  //
490  if ( lp_->getLeftPrec() != Teuchos::null ) {
491  lp_->applyLeftPrec( *R_, *Z_ );
492  if ( lp_->getRightPrec() != Teuchos::null ) {
493  Teuchos::RCP<MV> tmp = MVT::CloneCopy( *Z_);
494  lp_->applyRightPrec( *tmp, *Z_ );
495  }
496  }
497  else if ( lp_->getRightPrec() != Teuchos::null ) {
498  lp_->applyRightPrec( *R_, *Z_ );
499  }
500  else {
501  MVT::Assign( *R_, *Z_ );
502  }
503  //
504  if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() == Belos::TwoNorm) {
505  MVT::MvTransMv( one, *R_, *S_, rHs );
506  rHr_ = rHs(0,0);
507  rHz[0] = rHs(0,1);
508  } else
509  MVT::MvDot( *R_, *Z_, rHz );
510  //
511  beta[0] = rHz[0] / rHz_old[0];
512  //
513  MVT::MvAddMv( one, *Z_, beta[0], *P_, *P_ );
514 
515  // Condition estimate (if needed)
516  if (doCondEst_) {
517  if (iter_ > 1) {
518  diag_[iter_-1] = Teuchos::ScalarTraits<ScalarType>::real((beta_old * beta_old * pAp_old + pAp[0]) / rHz_old[0]);
519  offdiag_[iter_-2] = -Teuchos::ScalarTraits<ScalarType>::real(beta_old * pAp_old / (sqrt( rHz_old[0] * rHz_old2)));
520  }
521  else {
522  diag_[iter_-1] = Teuchos::ScalarTraits<ScalarType>::real(pAp[0] / rHz_old[0]);
523  }
524  rHz_old2 = rHz_old[0];
525  beta_old = beta[0];
526  pAp_old = pAp[0];
527  }
528 
529  } // end while (sTest_->checkStatus(this) != Passed)
530  }
531 
532 } // end Belos namespace
533 
534 #endif /* BELOS_CG_ITER_HPP */
Teuchos::RCP< const MV > R
The current residual.
Collection of types and exceptions used within the Belos solvers.
Teuchos::ScalarTraits< ScalarType > SCT
Definition: BelosCGIter.hpp:56
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
This class implements the preconditioned Conjugate Gradient (CG) iteration.
Definition: BelosCGIter.hpp:47
Structure to contain pointers to CGIteration state variables.
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
T & get(const std::string &name, T def_value)
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static magnitudeType real(T a)
Declaration of basic traits for the multivector type.
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
An implementation of StatusTestResNorm using a family of residual norms.
MultiVecTraits< ScalarType, MV > MVT
Definition: BelosCGIter.hpp:54
size_type size() const
A pure virtual class for defining the status tests for the Belos iterative solvers.
CGIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > &convTester, Teuchos::ParameterList &params)
CGIter constructor with linear problem, solver utilities, and parameter list of solver options...
Class which defines basic traits for the operator type.
CGIterateFailure is thrown when the CGIteration object is unable to compute the next iterate in the C...
Teuchos::RCP< const MV > AP
The matrix A applied to current decent direction vector.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
bool isInitialized()
States whether the solver has been initialized or not.
OperatorTraits< ScalarType, MV, OP > OPT
Definition: BelosCGIter.hpp:55
int getNumIters() const
Get the current iteration count.
Traits class which defines basic operations on multivectors.
CGPositiveDefiniteFailure is thrown when the the CG &#39;alpha = p^H*A*P&#39; value is less than zero...
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
void iterate()
This method performs CG iterations until the status test indicates the need to stop or an error occur...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation.
Teuchos::RCP< const MV > Z
The current preconditioned residual.
virtual ~CGIter()
Destructor.
Definition: BelosCGIter.hpp:74
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
SCT::magnitudeType MagnitudeType
Definition: BelosCGIter.hpp:57
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
void resetNumIters(int iter=0)
Reset the iteration count.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Class which defines basic traits for the operator type.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< const MV > P
The current decent direction vector.
void initializeCG(CGIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.

Generated on Fri Dec 20 2024 09:24:48 for Belos by doxygen 1.8.5