Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosPCPGIter.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_PCPG_ITER_HPP
11 #define BELOS_PCPG_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 
20 #include "BelosLinearProblem.hpp"
21 #include "BelosMatOrthoManager.hpp"
22 #include "BelosOutputManager.hpp"
23 #include "BelosStatusTest.hpp"
24 #include "BelosOperatorTraits.hpp"
25 #include "BelosMultiVecTraits.hpp"
26 #include "BelosCGIteration.hpp"
27 
30 #include "Teuchos_ScalarTraits.hpp"
32 #include "Teuchos_TimeMonitor.hpp"
33 
45 namespace Belos {
46 
48 
49 
54  template <class ScalarType, class MV>
55  struct PCPGIterState {
61  int curDim;
63  int prevUdim;
64 
67 
70 
73 
76 
79 
82 
88 
90  prevUdim(0),
91  R(Teuchos::null), Z(Teuchos::null),
92  P(Teuchos::null), AP(Teuchos::null),
93  U(Teuchos::null), C(Teuchos::null),
94  D(Teuchos::null)
95  {}
96  };
97 
99 
100  template<class ScalarType, class MV, class OP>
101  class PCPGIter : virtual public Iteration<ScalarType,MV,OP> {
102 
103  public:
104 
105  //
106  // Convenience typedefs
107  //
112 
114 
115 
124  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
127  Teuchos::ParameterList &params );
128 
130  virtual ~PCPGIter() {};
132 
133 
135 
136 
155  void iterate();
156 
176 
180  void initialize()
181  {
183  initialize(empty);
184  }
185 
195  state.Z = Z_; // CG state
196  state.P = P_;
197  state.AP = AP_;
198  state.R = R_;
199  state.U = U_; // seed state
200  state.C = C_;
201  state.D = D_;
202  state.curDim = curDim_;
203  state.prevUdim = prevUdim_;
204  return state;
205  }
206 
208 
209 
211 
212 
214  int getNumIters() const { return iter_; }
215 
217  void resetNumIters( int iter = 0 ) { iter_ = iter; }
218 
221  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> * /* norms */ ) const { return R_; }
222 
224 
227  Teuchos::RCP<MV> getCurrentUpdate() const { return Teuchos::null; }
228 
230  int getCurSubspaceDim() const {
231  if (!initialized_) return 0;
232  return curDim_;
233  };
234 
236  int getPrevSubspaceDim() const {
237  if (!initialized_) return 0;
238  return prevUdim_;
239  };
240 
242 
243 
245 
246 
248  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
249 
251  int getBlockSize() const { return 1; }
252 
254  int getNumRecycledBlocks() const { return savedBlocks_; }
255 
257 
259  void setBlockSize(int blockSize) {
260  TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
261  "Belos::PCPGIter::setBlockSize(): Cannot use a block size that is not one.");
262  }
263 
265  void setSize( int savedBlocks );
266 
268  bool isInitialized() { return initialized_; }
269 
271  void resetState();
272 
274 
275  private:
276 
277  //
278  // Internal methods
279  //
281  void setStateSize();
282 
283  //
284  // Classes inputed through constructor that define the linear problem to be solved.
285  //
290 
291  //
292  // Algorithmic parameters
293  // savedBlocks_ is the number of blocks allocated for the reused subspace
294  int savedBlocks_;
295  //
296  //
297  // Current solver state
298  //
299  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
300  // is capable of running; _initialize is controlled by the initialize() member method
301  // For the implications of the state of initialized_, please see documentation for initialize()
302  bool initialized_;
303 
304  // stateStorageInitialized_ indicates that the state storage has be initialized to the current
305  // savedBlocks_. State storage initialization may be postponed if the linear problem was
306  // generated without either the right-hand side or solution vectors.
307  bool stateStorageInitialized_;
308 
309  // keepDiagonal_ specifies that the iteration must keep the diagonal matrix of pivots
310  bool keepDiagonal_;
311 
312  // initDiagonal_ specifies that the iteration will reinitialize the diagonal matrix by zeroing
313  // out all entries before an iteration is started.
314  bool initDiagonal_;
315 
316  // Current subspace dimension
317  int curDim_;
318 
319  // Dimension of seed space used to solve current linear system
320  int prevUdim_;
321 
322  // Number of iterations performed
323  int iter_;
324  //
325  // State Storage ... of course this part is different for CG
326  //
327  // Residual
328  Teuchos::RCP<MV> R_;
329  //
330  // Preconditioned residual
331  Teuchos::RCP<MV> Z_;
332  //
333  // Direction std::vector
334  Teuchos::RCP<MV> P_;
335  //
336  // Operator applied to direction std::vector
337  Teuchos::RCP<MV> AP_;
338  //
339  // Recycled subspace vectors.
340  Teuchos::RCP<MV> U_;
341  //
342  // C = A * U, linear system is Ax=b
343  Teuchos::RCP<MV> C_;
344  //
345  // Projected matrices
346  // D_ : Diagonal matrix of pivots D = P'AP
348  };
349 
351  // Constructor.
352  template<class ScalarType, class MV, class OP>
354  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
357  Teuchos::ParameterList &params ):
358  lp_(problem),
359  om_(printer),
360  stest_(tester),
361  ortho_(ortho),
362  savedBlocks_(0),
363  initialized_(false),
364  stateStorageInitialized_(false),
365  keepDiagonal_(false),
366  initDiagonal_(false),
367  curDim_(0),
368  prevUdim_(0),
369  iter_(0)
370  {
371  // Get the maximum number of blocks allowed for this Krylov subspace
372 
373  TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter("Saved Blocks"), std::invalid_argument,
374  "Belos::PCPGIter::constructor: mandatory parameter \"Saved Blocks\" is not specified.");
375  int rb = Teuchos::getParameter<int>(params, "Saved Blocks");
376 
377  // Find out whether we are saving the Diagonal matrix.
378  keepDiagonal_ = params.get("Keep Diagonal", false);
379 
380  // Find out whether we are initializing the Diagonal matrix.
381  initDiagonal_ = params.get("Initialize Diagonal", false);
382 
383  // Set the number of blocks and allocate data
384  setSize( rb );
385  }
386 
388  // Set the block size and adjust as necessary
389  template <class ScalarType, class MV, class OP>
390  void PCPGIter<ScalarType,MV,OP>::setSize( int savedBlocks )
391  {
392  // allocate space only; perform no computation
393  // Any change in size invalidates the state of the solver as implemented here.
394 
395  TEUCHOS_TEST_FOR_EXCEPTION(savedBlocks <= 0, std::invalid_argument, "Belos::PCPGIter::setSize() was passed a non-positive argument for \"Num Saved Blocks\".");
396 
397  if ( savedBlocks_ != savedBlocks) {
398  stateStorageInitialized_ = false;
399  savedBlocks_ = savedBlocks;
400  initialized_ = false;
401  curDim_ = 0;
402  prevUdim_ = 0;
403  setStateSize(); // Use the current savedBlocks_ to initialize the state storage.
404  }
405  }
406 
408  // Enable the reuse of a single solver object for completely different linear systems
409  template <class ScalarType, class MV, class OP>
411  {
412  stateStorageInitialized_ = false;
413  initialized_ = false;
414  curDim_ = 0;
415  prevUdim_ = 0;
416  setStateSize();
417  }
418 
420  // Setup the state storage. Called by either initialize or, if savedBlocks_ changes, setSize.
421  template <class ScalarType, class MV, class OP>
423  {
424  if (!stateStorageInitialized_) {
425 
426  // Check if there is any multivector to clone from.
427  Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
428  Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
429  if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
430  return; // postpone exception
431  }
432  else {
433 
435  // blockSize*recycledBlocks dependent
436  int newsd = savedBlocks_ ; //int newsd = blockSize_* savedBlocks_ ;
437  //
438  // Initialize the CG state storage
439  // If the subspace is not initialized, generate it using the LHS or RHS from lp_.
440  // Generate CG state only if it does not exist, otherwise resize it.
441  if (Z_ == Teuchos::null) {
442  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
443  Z_ = MVT::Clone( *tmp, 1 );
444  }
445  if (P_ == Teuchos::null) {
446  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
447  P_ = MVT::Clone( *tmp, 1 );
448  }
449  if (AP_ == Teuchos::null) {
450  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
451  AP_ = MVT::Clone( *tmp, 1 );
452  }
453 
454  if (C_ == Teuchos::null) {
455 
456  // Get the multivector that is not null.
457  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
458  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
459  "Belos::PCPGIter::setStateSize(): linear problem does not specify multivectors to clone from.");
460  TEUCHOS_TEST_FOR_EXCEPTION( 0 != prevUdim_,std::invalid_argument,
461  "Belos::PCPGIter::setStateSize(): prevUdim not zero and C is null.");
462  C_ = MVT::Clone( *tmp, savedBlocks_ );
463  }
464  else {
465  // Generate C_ by cloning itself ONLY if more space is needed.
466  if (MVT::GetNumberVecs(*C_) < savedBlocks_ ) {
467  Teuchos::RCP<const MV> tmp = C_;
468  C_ = MVT::Clone( *tmp, savedBlocks_ );
469  }
470  }
471  if (U_ == Teuchos::null) {
472  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
473  TEUCHOS_TEST_FOR_EXCEPTION( 0 != prevUdim_,std::invalid_argument,
474  "Belos::PCPGIter::setStateSize(): prevUdim not zero and U is null.");
475  U_ = MVT::Clone( *tmp, savedBlocks_ );
476  }
477  else {
478  // Generate U_ by cloning itself ONLY if more space is needed.
479  if (MVT::GetNumberVecs(*U_) < savedBlocks_ ) {
480  Teuchos::RCP<const MV> tmp = U_;
481  U_ = MVT::Clone( *tmp, savedBlocks_ );
482  }
483  }
484  if (keepDiagonal_) {
485  if (D_ == Teuchos::null) {
487  }
488  if (initDiagonal_) {
489  D_->shape( newsd, newsd );
490  }
491  else {
492  if (D_->numRows() < newsd || D_->numCols() < newsd) {
493  D_->shapeUninitialized( newsd, newsd );
494  }
495  }
496  }
497  // State storage has now been initialized.
498  stateStorageInitialized_ = true;
499  } // if there is a vector to clone from
500  } // if !stateStorageInitialized_
501  } // end of setStateSize
502 
504  // Initialize the iteration object
505  template <class ScalarType, class MV, class OP>
507  {
508 
509  TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
510  "Belos::PCPGIter::initialize(): Cannot initialize state storage!");
511 
512  // Requirements: R_ and consistent multivectors widths and lengths
513  //
514  std::string errstr("Belos::PCPGIter::initialize(): Specified multivectors must have a consistent length and width.");
515 
516  if (newstate.R != Teuchos::null){
517 
518  R_ = newstate.R; // SolverManager::R_ == newstate.R == Iterator::R_
519  if (newstate.U == Teuchos::null){
520  prevUdim_ = 0;
521  newstate.U = U_;
522  newstate.C = C_;
523  }
524  else {
525  prevUdim_ = newstate.curDim;
526  if (newstate.C == Teuchos::null){ // Stub for new feature
527  std::vector<int> index(prevUdim_);
528  for (int i=0; i< prevUdim_; ++i)
529  index[i] = i;
530  Teuchos::RCP<const MV> Ukeff = MVT::CloneView( *newstate.U, index );
531  newstate.C = MVT::Clone( *newstate.U, prevUdim_ );
532  Teuchos::RCP<MV> Ckeff = MVT::CloneViewNonConst( *newstate.C, index );
533  lp_->apply( *Ukeff, *Ckeff );
534  }
535  curDim_ = prevUdim_ ;
536  }
537 
538  // Initialize the state storage if not already allocated in the constructor
539  if (!stateStorageInitialized_)
540  setStateSize();
541 
542  //TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_), std::invalid_argument, errstr );
543  //TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < 1, std::invalid_argument, errstr );
544 
545  newstate.prevUdim = prevUdim_; // big change in functionality from GCRODR
546  newstate.curDim = curDim_;
547 
548  //TEUCHOS_TEST_FOR_EXCEPTION(newstate.z->numRows() < curDim_ || newstate.z->numCols() < 1, std::invalid_argument, errstr);
549 
550  std::vector<int> zero_index(1);
551  zero_index[0] = 0;
552  if ( lp_->getLeftPrec() != Teuchos::null ) { // Compute the initial search direction
553  lp_->applyLeftPrec( *R_, *Z_ );
554  MVT::SetBlock( *Z_, zero_index , *P_ ); // P(:,zero_index) := Z
555  } else {
556  Z_ = R_;
557  MVT::SetBlock( *R_, zero_index, *P_ );
558  }
559 
560  std::vector<int> nextind(1);
561  nextind[0] = curDim_;
562 
563  MVT::SetBlock( *P_, nextind, *newstate.U ); // U(:,curDim_ ) := P_
564 
565  ++curDim_;
566  newstate.curDim = curDim_;
567 
568  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.U) != savedBlocks_ ,
569  std::invalid_argument, errstr );
570  if (newstate.U != U_) { // Why this is needed?
571  U_ = newstate.U;
572  }
573 
574  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.C) != savedBlocks_ ,
575  std::invalid_argument, errstr );
576  if (newstate.C != C_) {
577  C_ = newstate.C;
578  }
579  }
580  else {
581 
582  TEUCHOS_TEST_FOR_EXCEPTION(newstate.R == Teuchos::null,std::invalid_argument,
583  "Belos::PCPGIter::initialize(): PCPGStateIterState does not have initial kernel R_0.");
584  }
585 
586  // the solver is initialized
587  initialized_ = true;
588  }
589 
590 
592  // Iterate until the status test informs us we should stop.
593  template <class ScalarType, class MV, class OP>
595  {
596  //
597  // Allocate/initialize data structures
598  //
599  if (initialized_ == false) {
600  initialize();
601  }
602  const bool debug = false;
603 
604  // Allocate memory for scalars.
608  Teuchos::SerialDenseMatrix<int,ScalarType> rHz( 1, 1 ), rHz_old( 1, 1 );
609 
610  if( iter_ != 0 )
611  std::cout << " Iterate Warning: begin from nonzero iter_ ?" << std::endl; //DMD
612 
613  // GenOrtho Project Stubs
614  std::vector<int> prevInd;
618 
619  if( prevUdim_ ){
620  prevInd.resize( prevUdim_ );
621  for( int i=0; i<prevUdim_ ; i++) prevInd[i] = i;
622  CZ.reshape( prevUdim_ , 1 );
623  Uprev = MVT::CloneView(*U_, prevInd);
624  Cprev = MVT::CloneView(*C_, prevInd);
625  }
626 
627  // Get the current solution std::vector.
628  Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
629 
630  // Check that the current solution std::vector only has one column.
631  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*cur_soln_vec) != 1, CGIterationInitFailure,
632  "Belos::PCPGIter::iterate(): current linear system has more than one std::vector!" );
633 
634  //Check that the input is correctly set up
635  TEUCHOS_TEST_FOR_EXCEPTION( curDim_ != prevUdim_ + 1, CGIterationInitFailure,
636  "Belos::PCPGIter::iterate(): mistake in initialization !" );
637 
638 
639  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
640  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
641 
642 
643  std::vector<int> curind(1);
644  std::vector<ScalarType> rnorm(MVT::GetNumberVecs(*cur_soln_vec));
645  if (prevUdim_ > 0){ // A-orthonalize P=Z to Uprev
646  Teuchos::RCP<MV> P;
647  curind[0] = curDim_ - 1; // column = dimension - 1
648  P = MVT::CloneViewNonConst(*U_,curind);
649  MVT::MvTransMv( one, *Cprev, *P, CZ );
650  MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P ); // P -= U*(C'Z)
651 
652  if( debug ){
653  MVT::MvTransMv( one, *Cprev, *P, CZ );
654  std::cout << " Input CZ post ortho " << std::endl;
655  CZ.print( std::cout );
656  }
657  if( curDim_ == savedBlocks_ ){
658  std::vector<int> zero_index(1);
659  zero_index[0] = 0;
660  MVT::SetBlock( *P, zero_index, *P_ );
661  }
662  P = Teuchos::null;
663  }
664 
665  // Compute first <r,z> a.k.a. rHz
666  MVT::MvTransMv( one, *R_, *Z_, rHz );
667 
669  // iterate until the status test is satisfied
670  //
671  while (stest_->checkStatus(this) != Passed ) {
673  Teuchos::RCP<MV> AP;
674  iter_++; // The next iteration begins.
675  //std::vector<int> curind(1);
676  curind[0] = curDim_ - 1; // column = dimension - 1
677  if( debug ){
678  MVT::MvNorm(*R_, rnorm);
679  std::cout << iter_ << " " << curDim_ << " " << rnorm[0] << std::endl;
680  }
681  if( prevUdim_ + iter_ < savedBlocks_ ){
682  P = MVT::CloneView(*U_,curind);
683  AP = MVT::CloneViewNonConst(*C_,curind);
684  lp_->applyOp( *P, *AP );
685  MVT::MvTransMv( one, *P, *AP, pAp );
686  }else{
687  if( prevUdim_ + iter_ == savedBlocks_ ){
688  AP = MVT::CloneViewNonConst(*C_,curind);
689  lp_->applyOp( *P_, *AP );
690  MVT::MvTransMv( one, *P_, *AP, pAp );
691  }else{
692  lp_->applyOp( *P_, *AP_ );
693  MVT::MvTransMv( one, *P_, *AP_, pAp );
694  }
695  }
696 
697  if( keepDiagonal_ && prevUdim_ + iter_ <= savedBlocks_ )
698  (*D_)(iter_ -1 ,iter_ -1 ) = pAp(0,0);
699 
700  // positive pAp required
702  "Belos::PCPGIter::iterate(): non-positive value for p^H*A*p encountered!" );
703 
704  // alpha := <R_,Z_> / <P,AP>
705  alpha(0,0) = rHz(0,0) / pAp(0,0);
706 
707  // positive alpha required
709  "Belos::PCPGIter::iterate(): non-positive value for alpha encountered!" );
710 
711  // solution update x += alpha * P
712  if( curDim_ < savedBlocks_ ){
713  MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P, *cur_soln_vec );
714  }else{
715  MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P_, *cur_soln_vec );
716  }
717  //lp_->updateSolution(); ... does nothing.
718  //
719  // The denominator of beta is saved before residual is updated [ old <R_, Z_> ].
720  //
721  rHz_old(0,0) = rHz(0,0);
722  //
723  // residual update R_ := R_ - alpha * AP
724  //
725  if( prevUdim_ + iter_ <= savedBlocks_ ){
726  MVT::MvAddMv( one, *R_, -alpha(0,0), *AP, *R_ );
727  AP = Teuchos::null;
728  }else{
729  MVT::MvAddMv( one, *R_, -alpha(0,0), *AP_, *R_ );
730  }
731  //
732  // update beta := [ new <R_, Z_> ] / [ old <R_, Z_> ] and the search direction p.
733  //
734  if ( lp_->getLeftPrec() != Teuchos::null ) {
735  lp_->applyLeftPrec( *R_, *Z_ );
736  } else {
737  Z_ = R_;
738  }
739  //
740  MVT::MvTransMv( one, *R_, *Z_, rHz );
741  //
742  beta(0,0) = rHz(0,0) / rHz_old(0,0);
743  //
744  if( curDim_ < savedBlocks_ ){
745  curDim_++; // update basis dim
746  curind[0] = curDim_ - 1;
747  Teuchos::RCP<MV> Pnext = MVT::CloneViewNonConst(*U_,curind);
748  MVT::MvAddMv( one, *Z_, beta(0,0), *P, *Pnext );
749  if( prevUdim_ ){ // Deflate seed space
750  MVT::MvTransMv( one, *Cprev, *Z_, CZ );
751  MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *Pnext ); // Pnext -= U*(C'Z)
752  if( debug ){
753  std::cout << " Check CZ " << std::endl;
754  MVT::MvTransMv( one, *Cprev, *Pnext, CZ );
755  CZ.print( std::cout );
756  }
757  }
758  P = Teuchos::null;
759  if( curDim_ == savedBlocks_ ){
760  std::vector<int> zero_index(1);
761  zero_index[0] = 0;
762  MVT::SetBlock( *Pnext, zero_index, *P_ );
763  }
764  Pnext = Teuchos::null;
765  }else{
766  MVT::MvAddMv( one, *Z_, beta(0,0), *P_, *P_ );
767  if( prevUdim_ ){ // Deflate seed space
768  MVT::MvTransMv( one, *Cprev, *Z_, CZ );
769  MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P_ ); // P_ -= U*(C'Z)
770 
771  if( debug ){
772  std::cout << " Check CZ " << std::endl;
773  MVT::MvTransMv( one, *Cprev, *P_, CZ );
774  CZ.print( std::cout );
775  }
776  }
777  }
778  // CGB: 5/26/2010
779  // this RCP<const MV> P was previously a variable outside the loop. however, it didn't appear to be see any use between
780  // loop iterations. therefore, I moved it inside to avoid scoping errors with previously used variables named P.
781  // to ensure that this wasn't a bug, I verify below that we have set P == null, i.e., that we are not going to use it again
782  // same for AP
783  TEUCHOS_TEST_FOR_EXCEPTION( AP != Teuchos::null || P != Teuchos::null, std::logic_error, "Loop recurrence violated. Please contact Belos team.");
784  } // end coupled two-term recursion
785  if( prevUdim_ + iter_ < savedBlocks_ ) --curDim_; // discard negligible search direction
786  }
787 
788 } // end Belos namespace
789 
790 #endif /* BELOS_PCPG_ITER_HPP */
Collection of types and exceptions used within the Belos solvers.
int prevUdim
Number of block columns in matrices C and U before current iteration.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
bool isInitialized()
States whether the solver has been initialized or not.
Class which manages the output and verbosity of the Belos solvers.
Teuchos::RCP< MV > R
The current residual.
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
PCPGIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Pure virtual base class for defining the status testing capabilities of Belos.
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
int getCurSubspaceDim() const
Get the current dimension of the whole seed subspace.
Declaration of basic traits for the multivector type.
virtual std::ostream & print(std::ostream &os) const
void setSize(int savedBlocks)
Set the maximum number of saved or recycled blocks used by the iterative solver.
int curDim
The current dimension of the reduction.
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::RCP< MV > P
The current decent direction std::vector.
Teuchos::RCP< MV > AP
The matrix A applied to current decent direction std::vector.
Structure to contain pointers to PCPGIter state variables.
OperatorTraits< ScalarType, MV, OP > OPT
int getPrevSubspaceDim() const
Get the dimension of the search subspace used to solve the current solution to the linear problem...
A pure virtual class for defining the status tests for the Belos iterative solvers.
int getBlockSize() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
Class which defines basic traits for the operator type.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.
CGPositiveDefiniteFailure is thrown when the the CG &#39;alpha = p^H*A*P&#39; value is less than zero...
bool isParameter(const std::string &name) const
virtual ~PCPGIter()
Destructor.
void setBlockSize(int blockSize)
Get the blocksize to be used by the iterative solver in solving this linear problem.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int getNumIters() const
Get the current iteration count.
void initialize()
Initialize the solver with the initial vectors from the linear problem. An exception is thrown if ini...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getNumRecycledBlocks() const
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem...
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
SCT::magnitudeType MagnitudeType
Teuchos::RCP< MV > U
The recycled subspace.
void resetState()
tell the Iterator to &quot;reset&quot; itself; delete and rebuild the seed space.
int reshape(OrdinalType numRows, OrdinalType numCols)
void iterate()
PCPGIter iterates CG until the status test either requests a stop or detects an error. In the latter case, std::exception is thrown.
Teuchos::RCP< MV > Z
The current preconditioned residual.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
Get the norms of the residuals native to the solver.
Class which defines basic traits for the operator type.
CGIterationInitFailure is thrown when the CGIteration object is unable to generate an initial iterate...
PCPGIter(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< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
PCPGIter constructor with linear problem, solver utilities, and parameter list of solver options...
Belos&#39;s templated virtual class for providing routines for orthogonalization and orthonormzalition of...
This class implements the PCPG iteration, where a single-std::vector Krylov subspace is constructed...
Belos header file which uses auto-configuration information to include necessary C++ headers...
MultiVecTraits< ScalarType, MV > MVT
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system solution?.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Teuchos::ScalarTraits< ScalarType > SCT

Generated on Fri Nov 22 2024 09:23:06 for Belos by doxygen 1.8.5