Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosPseudoBlockCGIter.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_PSEUDO_BLOCK_CG_ITER_HPP
11 #define BELOS_PSEUDO_BLOCK_CG_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosCGIteration.hpp"
20 
21 #include "BelosLinearProblem.hpp"
22 #include "BelosMatOrthoManager.hpp"
23 #include "BelosOutputManager.hpp"
24 #include "BelosStatusTest.hpp"
25 #include "BelosOperatorTraits.hpp"
26 #include "BelosMultiVecTraits.hpp"
27 
30 #include "Teuchos_ScalarTraits.hpp"
32 #include "Teuchos_TimeMonitor.hpp"
33 
45 namespace Belos {
46 
47  template<class ScalarType, class MV, class OP>
48  class PseudoBlockCGIter : virtual public CGIteration<ScalarType,MV,OP> {
49 
50  public:
51 
52  //
53  // Convenience typedefs
54  //
59 
61 
62 
69  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
71  Teuchos::ParameterList &params );
72 
74  virtual ~PseudoBlockCGIter() {};
76 
77 
79 
80 
94  void iterate();
95 
117 
121  void initialize()
122  {
124  initializeCG(empty);
125  }
126 
136  state.R = R_;
137  state.P = P_;
138  state.AP = AP_;
139  state.Z = Z_;
140  return state;
141  }
142 
144 
145 
147 
148 
150  int getNumIters() const { return iter_; }
151 
153  void resetNumIters( int iter = 0 ) { iter_ = iter; }
154 
157  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> * /* norms */ ) const { return R_; }
158 
160 
162  Teuchos::RCP<MV> getCurrentUpdate() const { return Teuchos::null; }
163 
165 
167 
168 
170  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
171 
173  int getBlockSize() const { return 1; }
174 
176  void setBlockSize(int blockSize) {
177  TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
178  "Belos::PseudoBlockCGIter::setBlockSize(): Cannot use a block size that is not one.");
179  }
180 
182  bool isInitialized() { return initialized_; }
183 
185 
187  void setDoCondEst(bool val) {
188  if (numEntriesForCondEst_) doCondEst_=val;
189  }
190 
193  // NOTE (mfh 30 Jul 2015) See note on getOffDiag() below.
194  // getDiag() didn't actually throw for me in that case, but why
195  // not be cautious?
196  typedef typename Teuchos::ArrayView<MagnitudeType>::size_type size_type;
197  if (static_cast<size_type> (iter_) >= diag_.size ()) {
198  return diag_ ();
199  } else {
200  return diag_ (0, iter_);
201  }
202  }
203 
206  // NOTE (mfh 30 Jul 2015) The implementation as I found it
207  // returned "offdiag(0,iter_)". This breaks (Teuchos throws in
208  // debug mode) when the maximum number of iterations has been
209  // reached, because iter_ == offdiag_.size() in that case. The
210  // new logic fixes this.
211  typedef typename Teuchos::ArrayView<MagnitudeType>::size_type size_type;
212  if (static_cast<size_type> (iter_) >= offdiag_.size ()) {
213  return offdiag_ ();
214  } else {
215  return offdiag_ (0, iter_);
216  }
217  }
218 
219  private:
220 
221  //
222  // Classes inputed through constructor that define the linear problem to be solved.
223  //
227 
228  //
229  // Algorithmic parameters
230  //
231  // numRHS_ is the current number of linear systems being solved.
232  int numRHS_;
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  // Current number of iterations performed.
243  int iter_;
244 
245  // Assert that the matrix is positive definite
246  bool assertPositiveDefiniteness_;
247 
248  // Tridiagonal system for condition estimation (if needed)
249  Teuchos::ArrayRCP<MagnitudeType> diag_, offdiag_;
250  ScalarType pAp_old_, beta_old_, rHz_old2_; // Put scalars here so that estimate is correct for multiple RHS, when deflation occurs.
251  int numEntriesForCondEst_;
252  bool doCondEst_;
253 
254  //
255  // State Storage
256  //
257  // Residual
258  Teuchos::RCP<MV> R_;
259  //
260  // Preconditioned residual
261  Teuchos::RCP<MV> Z_;
262  //
263  // Direction vector
264  Teuchos::RCP<MV> P_;
265  //
266  // Operator applied to direction vector
267  Teuchos::RCP<MV> AP_;
268 
269  };
270 
272  // Constructor.
273  template<class ScalarType, class MV, class OP>
275  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
277  Teuchos::ParameterList &params ):
278  lp_(problem),
279  om_(printer),
280  stest_(tester),
281  numRHS_(0),
282  initialized_(false),
283  iter_(0),
284  assertPositiveDefiniteness_( params.get("Assert Positive Definiteness", true) ),
285  numEntriesForCondEst_(params.get("Max Size For Condest",0) ),
286  doCondEst_(false)
287  {
288  }
289 
290 
292  // Initialize this iteration object
293  template <class ScalarType, class MV, class OP>
295  {
296  // Check if there is any mltivector to clone from.
297  Teuchos::RCP<const MV> lhsMV = lp_->getCurrLHSVec();
298  Teuchos::RCP<const MV> rhsMV = lp_->getCurrRHSVec();
299  TEUCHOS_TEST_FOR_EXCEPTION((lhsMV==Teuchos::null && rhsMV==Teuchos::null),std::invalid_argument,
300  "Belos::PseudoBlockCGIter::initialize(): Cannot initialize state storage!");
301 
302  // Get the multivector that is not null.
303  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
304 
305  // Get the number of right-hand sides we're solving for now.
306  int numRHS = MVT::GetNumberVecs(*tmp);
307  numRHS_ = numRHS;
308 
309  // Initialize the state storage
310  // If the subspace has not be initialized before or has changed sizes, generate it using the LHS or RHS from lp_.
311  if (Teuchos::is_null(R_) || MVT::GetNumberVecs(*R_)!=numRHS_) {
312  R_ = MVT::Clone( *tmp, numRHS_ );
313  Z_ = MVT::Clone( *tmp, numRHS_ );
314  P_ = MVT::Clone( *tmp, numRHS_ );
315  AP_ = MVT::Clone( *tmp, numRHS_ );
316  }
317 
318  // Tracking information for condition number estimation
319  if(numEntriesForCondEst_ > 0) {
320  diag_.resize(numEntriesForCondEst_);
321  offdiag_.resize(numEntriesForCondEst_-1);
322  }
323 
324  // NOTE: In CGIter R_, the initial residual, is required!!!
325  //
326  std::string errstr("Belos::BlockPseudoCGIter::initialize(): Specified multivectors must have a consistent length and width.");
327 
328  // Create convenience variables for zero and one.
329  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
331 
332  if (!Teuchos::is_null(newstate.R)) {
333 
334  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*R_),
335  std::invalid_argument, errstr );
336  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.R) != numRHS_,
337  std::invalid_argument, errstr );
338 
339  // Copy basis vectors from newstate into V
340  if (newstate.R != R_) {
341  // copy over the initial residual (unpreconditioned).
342  MVT::MvAddMv( one, *newstate.R, zero, *newstate.R, *R_ );
343  }
344 
345  // Compute initial direction vectors
346  // Initially, they are set to the preconditioned residuals
347  //
348  if ( lp_->getLeftPrec() != Teuchos::null ) {
349  lp_->applyLeftPrec( *R_, *Z_ );
350  if ( lp_->getRightPrec() != Teuchos::null ) {
351  Teuchos::RCP<MV> tmp1 = MVT::Clone( *Z_, numRHS_ );
352  lp_->applyRightPrec( *Z_, *tmp1 );
353  Z_ = tmp1;
354  }
355  }
356  else if ( lp_->getRightPrec() != Teuchos::null ) {
357  lp_->applyRightPrec( *R_, *Z_ );
358  }
359  else {
360  Z_ = R_;
361  }
362  MVT::MvAddMv( one, *Z_, zero, *Z_, *P_ );
363  }
364  else {
365 
366  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::is_null(newstate.R),std::invalid_argument,
367  "Belos::CGIter::initialize(): CGStateIterState does not have initial residual.");
368  }
369 
370  // The solver is initialized
371  initialized_ = true;
372  }
373 
374 
376  // Iterate until the status test informs us we should stop.
377  template <class ScalarType, class MV, class OP>
379  {
380  //
381  // Allocate/initialize data structures
382  //
383  if (initialized_ == false) {
384  initialize();
385  }
386 
387  // Allocate memory for scalars.
388  int i=0;
389  std::vector<int> index(1);
390  std::vector<ScalarType> rHz( numRHS_ ), rHz_old( numRHS_ ), pAp( numRHS_ ), beta( numRHS_ );
391  Teuchos::SerialDenseMatrix<int, ScalarType> alpha( numRHS_,numRHS_ );
392 
393  // Create convenience variables for zero and one.
394  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
396 
397  // Get the current solution std::vector.
398  Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
399 
400  // Compute first <r,z> a.k.a. rHz
401  MVT::MvDot( *R_, *Z_, rHz );
402 
403  if ( assertPositiveDefiniteness_ )
404  for (i=0; i<numRHS_; ++i)
405  TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(rHz[i]) < zero,
407  "Belos::PseudoBlockCGIter::iterate(): negative value for r^H*M*r encountered!" );
408 
410  // Iterate until the status test tells us to stop.
411  //
412  while (stest_->checkStatus(this) != Passed) {
413 
414  // Increment the iteration
415  iter_++;
416 
417  // Multiply the current direction std::vector by A and store in AP_
418  lp_->applyOp( *P_, *AP_ );
419 
420  // Compute alpha := <R_,Z_> / <P_,AP_>
421  MVT::MvDot( *P_, *AP_, pAp );
422 
423  for (i=0; i<numRHS_; ++i) {
424  if ( assertPositiveDefiniteness_ )
425  // Check that pAp[i] is a positive number!
426  TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(pAp[i]) <= zero,
428  "Belos::PseudoBlockCGIter::iterate(): non-positive value for p^H*A*p encountered!" );
429 
430  alpha(i,i) = rHz[i] / pAp[i];
431  }
432 
433  //
434  // Update the solution std::vector x := x + alpha * P_
435  //
436  MVT::MvTimesMatAddMv( one, *P_, alpha, one, *cur_soln_vec );
437  lp_->updateSolution();// what does this do?
438  //
439  // Save the denominator of beta before residual is updated [ old <R_, Z_> ]
440  //
441  for (i=0; i<numRHS_; ++i) {
442  rHz_old[i] = rHz[i];
443  }
444  //
445  // Compute the new residual R_ := R_ - alpha * AP_
446  //
447  MVT::MvTimesMatAddMv( -one, *AP_, alpha, one, *R_ );
448  //
449  // Compute beta := [ new <R_, Z_> ] / [ old <R_, Z_> ],
450  // and the new direction std::vector p.
451  //
452  if ( lp_->getLeftPrec() != Teuchos::null ) {
453  lp_->applyLeftPrec( *R_, *Z_ );
454  if ( lp_->getRightPrec() != Teuchos::null ) {
455  Teuchos::RCP<MV> tmp = MVT::Clone( *Z_, numRHS_ );
456  lp_->applyRightPrec( *Z_, *tmp );
457  Z_ = tmp;
458  }
459  }
460  else if ( lp_->getRightPrec() != Teuchos::null ) {
461  lp_->applyRightPrec( *R_, *Z_ );
462  }
463  else {
464  Z_ = R_;
465  }
466  //
467  MVT::MvDot( *R_, *Z_, rHz );
468  if ( assertPositiveDefiniteness_ )
469  for (i=0; i<numRHS_; ++i)
470  TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(rHz[i]) < zero,
472  "Belos::PseudoBlockCGIter::iterate(): negative value for r^H*M*r encountered!" );
473  //
474  // Update the search directions.
475  for (i=0; i<numRHS_; ++i) {
476  beta[i] = rHz[i] / rHz_old[i];
477  index[0] = i;
478  Teuchos::RCP<const MV> Z_i = MVT::CloneView( *Z_, index );
479  Teuchos::RCP<MV> P_i = MVT::CloneViewNonConst( *P_, index );
480  MVT::MvAddMv( one, *Z_i, beta[i], *P_i, *P_i );
481  }
482 
483  // Condition estimate (if needed)
484  if (doCondEst_ && (iter_ - 1) < diag_.size()) {
485  if (iter_ > 1) {
486  diag_[iter_-1] = Teuchos::ScalarTraits<ScalarType>::real((beta_old_ * beta_old_ * pAp_old_ + pAp[0]) / rHz_old[0]);
487  offdiag_[iter_-2] = -Teuchos::ScalarTraits<ScalarType>::real(beta_old_ * pAp_old_ / (sqrt( rHz_old[0] * rHz_old2_)));
488  }
489  else {
490  diag_[iter_-1] = Teuchos::ScalarTraits<ScalarType>::real(pAp[0] / rHz_old[0]);
491  }
492  rHz_old2_ = rHz_old[0];
493  beta_old_ = beta[0];
494  pAp_old_ = pAp[0];
495  }
496 
497 
498  //
499  } // end while (sTest_->checkStatus(this) != Passed)
500  }
501 
502 } // end Belos namespace
503 
504 #endif /* BELOS_PSEUDO_BLOCK_CG_ITER_HPP */
Teuchos::RCP< const MV > R
The current residual.
Collection of types and exceptions used within the Belos solvers.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
Get the norms of the residuals native to the solver.
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.
virtual ~PseudoBlockCGIter()
Destructor.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation.
Structure to contain pointers to CGIteration state variables.
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
bool is_null(const std::shared_ptr< T > &p)
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.
size_type size() const
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
Teuchos::RCP< const MV > AP
The matrix A applied to current decent direction vector.
OperatorTraits< ScalarType, MV, OP > OPT
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
CGPositiveDefiniteFailure is thrown when the the CG &#39;alpha = p^H*A*P&#39; value is less than zero...
void resetNumIters(int iter=0)
Reset the iteration count.
This class implements the pseudo-block CG iteration, where the basic CG algorithm is performed on all...
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
MultiVecTraits< ScalarType, MV > MVT
void iterate()
This method performs CG iterations on each linear system until the status test indicates the need to ...
A linear system to solve, and its associated information.
int getNumIters() const
Get the current iteration count.
Class which describes the linear problem to be solved by the iterative solver.
Teuchos::ScalarTraits< ScalarType > SCT
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Teuchos::RCP< const MV > Z
The current preconditioned residual.
PseudoBlockCGIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, Teuchos::ParameterList &params)
PseudoBlockCGIter constructor with linear problem, solver utilities, and parameter list of solver opt...
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation.
bool isInitialized()
States whether the solver has been initialized or not.
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.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...

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