Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosGCRODRIter.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_GCRODR_ITER_HPP
11 #define BELOS_GCRODR_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 
27 #include "Teuchos_BLAS.hpp"
30 #include "Teuchos_ScalarTraits.hpp"
32 #include "Teuchos_TimeMonitor.hpp"
33 
46 namespace Belos {
47 
49 
50 
55  template <class ScalarType, class MV>
56  struct GCRODRIterState {
61  int curDim;
62 
65 
68 
75 
79 
80  GCRODRIterState() : curDim(0), V(Teuchos::null),
81  U(Teuchos::null), C(Teuchos::null),
82  H(Teuchos::null), B(Teuchos::null)
83  {}
84  };
85 
87 
89 
90 
103  public:
104  GCRODRIterInitFailure(const std::string& what_arg) : BelosError(what_arg) {}
105  };
106 
114  public:
115  GCRODRIterOrthoFailure(const std::string& what_arg) : BelosError(what_arg) {}
116  };
117 
119 
120 
121  template<class ScalarType, class MV, class OP>
122  class GCRODRIter : virtual public Iteration<ScalarType,MV,OP> {
123 
124  public:
125 
126  //
127  // Convenience typedefs
128  //
133 
135 
136 
146  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
149  Teuchos::ParameterList &params );
150 
152  virtual ~GCRODRIter() {};
154 
155 
157 
158 
180  void iterate();
181 
204 
208  void initialize() {
210  initialize(empty);
211  }
212 
222  state.curDim = curDim_;
223  state.V = V_;
224  state.U = U_;
225  state.C = C_;
226  state.H = H_;
227  state.B = B_;
228  return state;
229  }
230 
232 
233 
235 
236 
238  int getNumIters() const { return iter_; }
239 
241  void resetNumIters( int iter = 0 ) { iter_ = iter; }
242 
245  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
246 
248 
254 
256 
259  void updateLSQR( int dim = -1 );
260 
262  int getCurSubspaceDim() const {
263  if (!initialized_) return 0;
264  return curDim_;
265  };
266 
268  int getMaxSubspaceDim() const { return numBlocks_; }
269 
271 
272 
274 
275 
277  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
278 
280  int getNumBlocks() const { return numBlocks_; }
281 
283  void setNumBlocks(int numBlocks) { setSize( recycledBlocks_, numBlocks ); };
284 
286  int getRecycledBlocks() const { return recycledBlocks_; }
287 
289  void setRecycledBlocks(int recycledBlocks) { setSize( recycledBlocks, numBlocks_ ); };
290 
292  int getBlockSize() const { return 1; }
293 
295  void setBlockSize(int blockSize) {
296  TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,"Belos::GCRODRIter::setBlockSize(): Cannot use a block size that is not one.");
297  }
298 
300  void setSize( int recycledBlocks, int numBlocks ) {
301  // only call resize if size changed
302  if ( (recycledBlocks_ != recycledBlocks) || (numBlocks_ != numBlocks) ) {
303  recycledBlocks_ = recycledBlocks;
304  numBlocks_ = numBlocks;
305  cs_.sizeUninitialized( numBlocks_+1 );
306  sn_.sizeUninitialized( numBlocks_+1 );
307  z_.sizeUninitialized( numBlocks_+1 );
308  R_.shapeUninitialized( numBlocks_+1,numBlocks );
309  }
310  }
311 
313  bool isInitialized() { return initialized_; }
314 
316 
317  private:
318 
319  //
320  // Internal methods
321  //
322 
323  //
324  // Classes inputed through constructor that define the linear problem to be solved.
325  //
330 
331  //
332  // Algorithmic parameters
333  //
334  // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
335  int numBlocks_;
336 
337  // recycledBlocks_ is the size of the allocated space for the recycled subspace, in blocks.
338  int recycledBlocks_;
339 
340  // Storage for QR factorization of the least squares system.
343 
344  //
345  // Current solver state
346  //
347  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
348  // is capable of running; _initialize is controlled by the initialize() member method
349  // For the implications of the state of initialized_, please see documentation for initialize()
350  bool initialized_;
351 
352  // Current subspace dimension, and number of iterations performed.
353  int curDim_, iter_;
354 
355  //
356  // State Storage
357  //
358  // Krylov vectors.
359  Teuchos::RCP<MV> V_;
360  //
361  // Recycled subspace vectors.
362  Teuchos::RCP<MV> U_, C_;
363  //
364  // Projected matrices
365  // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
367  //
368  // B_ : Projected matrix from the recycled subspace B = C^H*A*V
370  //
371  // QR decomposition of Projected matrices for solving the least squares system HY = B.
372  // R_: Upper triangular reduction of H
373  // z_: Q applied to right-hand side of the least squares system
376  };
377 
379  // Constructor.
380  template<class ScalarType, class MV, class OP>
382  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
385  Teuchos::ParameterList &params ):
386  lp_(problem), om_(printer), stest_(tester), ortho_(ortho) {
387  numBlocks_ = 0;
388  recycledBlocks_ = 0;
389  initialized_ = false;
390  curDim_ = 0;
391  iter_ = 0;
392  V_ = Teuchos::null;
393  U_ = Teuchos::null;
394  C_ = Teuchos::null;
395  H_ = Teuchos::null;
396  B_ = Teuchos::null;
397 
398  // Get the maximum number of blocks allowed for this Krylov subspace
399  TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter("Num Blocks"), std::invalid_argument, "Belos::GCRODRIter::constructor: mandatory parameter \"Num Blocks\" is not specified.");
400  int nb = Teuchos::getParameter<int>(params, "Num Blocks");
401 
402  TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter("Recycled Blocks"), std::invalid_argument,"Belos::GCRODRIter::constructor: mandatory parameter \"Recycled Blocks\" is not specified.");
403  int rb = Teuchos::getParameter<int>(params, "Recycled Blocks");
404 
405  TEUCHOS_TEST_FOR_EXCEPTION(nb <= 0, std::invalid_argument, "Belos::GCRODRIter() was passed a non-positive argument for \"Num Blocks\".");
406  TEUCHOS_TEST_FOR_EXCEPTION(rb >= nb, std::invalid_argument, "Belos::GCRODRIter() the number of recycled blocks is larger than the allowable subspace.");
407 
408  numBlocks_ = nb;
409  recycledBlocks_ = rb;
410 
411  cs_.sizeUninitialized( numBlocks_+1 );
412  sn_.sizeUninitialized( numBlocks_+1 );
413  z_.sizeUninitialized( numBlocks_+1 );
414  R_.shapeUninitialized( numBlocks_+1,numBlocks_ );
415 
416  }
417 
419  // Get the current update from this subspace.
420  template <class ScalarType, class MV, class OP>
422  //
423  // If this is the first iteration of the Arnoldi factorization,
424  // there is no update, so return Teuchos::null.
425  //
426  Teuchos::RCP<MV> currentUpdate = Teuchos::null;
427  if (curDim_==0) {
428  return currentUpdate;
429  } else {
430  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
431  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
433  currentUpdate = MVT::Clone( *V_, 1 );
434  //
435  // Make a view and then copy the RHS of the least squares problem. DON'T OVERWRITE IT!
436  //
438  //
439  // Solve the least squares problem.
440  //
442  Teuchos::NON_UNIT_DIAG, curDim_, 1, one,
443  R_.values(), R_.stride(), y.values(), y.stride() );
444  //
445  // Compute the current update from the Krylov basis; V(:,1:curDim_)*y.
446  //
447  std::vector<int> index(curDim_);
448  for ( int i=0; i<curDim_; i++ ) index[i] = i;
449  Teuchos::RCP<const MV> Vjp1 = MVT::CloneView( *V_, index );
450  MVT::MvTimesMatAddMv( one, *Vjp1, y, zero, *currentUpdate );
451  //
452  // Add in portion of update from recycled subspace U; U(:,1:recycledBlocks_)*B*y.
453  //
454  if (U_ != Teuchos::null) {
455  Teuchos::SerialDenseMatrix<int,ScalarType> z(recycledBlocks_,1);
456  Teuchos::SerialDenseMatrix<int,ScalarType> subB( Teuchos::View, *B_, recycledBlocks_, curDim_ );
457  z.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, subB, y, zero );
458  MVT::MvTimesMatAddMv( -one, *U_, z, one, *currentUpdate );
459  }
460  }
461  return currentUpdate;
462  }
463 
464 
466  // Get the native residuals stored in this iteration.
467  // Note: This method does not return a MultiVector of the residual vectors, only return Teuchos::null
468  template <class ScalarType, class MV, class OP>
470  //
471  // NOTE: Make sure the incoming std::vector is the correct size!
472  //
473  if ( norms && (int)norms->size()==0 )
474  norms->resize( 1 );
475 
476  if (norms) {
478  (*norms)[0] = blas.NRM2( 1, &z_(curDim_), 1);
479  }
480  return Teuchos::null;
481  }
482 
483 
484 
486  // Initialize this iteration object
487  template <class ScalarType, class MV, class OP>
489 
490  if (newstate.V != Teuchos::null && newstate.H != Teuchos::null) {
491  curDim_ = newstate.curDim;
492  V_ = newstate.V;
493  U_ = newstate.U;
494  C_ = newstate.C;
495  H_ = newstate.H;
496  B_ = newstate.B;
497  }
498  else {
499  TEUCHOS_TEST_FOR_EXCEPTION(newstate.V == Teuchos::null,std::invalid_argument,"Belos::GCRODRIter::initialize(): GCRODRIterState does not have V initialized.");
500  TEUCHOS_TEST_FOR_EXCEPTION(newstate.H == Teuchos::null,std::invalid_argument,"Belos::GCRODRIter::initialize(): GCRODRIterState does not have H initialized.");
501  }
502 
503  // the solver is initialized
504  initialized_ = true;
505 
506  }
507 
508 
510  // Iterate until the status test informs us we should stop.
511  template <class ScalarType, class MV, class OP>
513 
514  TEUCHOS_TEST_FOR_EXCEPTION( initialized_ == false, GCRODRIterInitFailure,"Belos::GCRODRIter::iterate(): GCRODRIter class not initialized." );
515 
516  // Force call to setsize to ensure internal storage is correct dimension
517  setSize( recycledBlocks_, numBlocks_ );
518 
519  Teuchos::RCP<MV> Vnext;
521  std::vector<int> curind(1);
522 
523  // z_ must be zeroed out in order to compute Givens rotations correctly
524  z_.putScalar(0.0);
525 
526  // Orthonormalize the new V_0
527  curind[0] = 0;
528  Vnext = MVT::CloneViewNonConst(*V_,curind);
529  // Orthonormalize first column
530  // First, get a monoelemental matrix to hold the orthonormalization coefficients
533  int rank = ortho_->normalize( *Vnext, z0 );
534  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1,GCRODRIterOrthoFailure, "Belos::GCRODRIter::iterate(): couldn't generate basis of full rank.");
535  // Copy the scalar coefficient back into the z_ vector
536  z_(0) = (*z0)(0,0);
537 
538  std::vector<int> prevind(numBlocks_+1);
539 
541  // iterate until the status test tells us to stop.
542  //
543  // also break if our basis is full
544  //
545  if (U_ == Teuchos::null) { // iterate without recycle space (because none is available yet)
546  while (stest_->checkStatus(this) != Passed && curDim_+1 <= numBlocks_) {
547  iter_++;
548  int lclDim = curDim_ + 1;
549 
550  // Get next index in basis
551  curind[0] = lclDim;
552  Vnext = MVT::CloneViewNonConst(*V_,curind);
553 
554  // Get previous index in basis
555  curind[0] = curDim_;
556  Vprev = MVT::CloneView(*V_,curind);
557 
558  // Compute the next vector in the Krylov basis: Vnext = Op*Vprev
559  lp_->apply(*Vprev,*Vnext);
560 
561  // Remove all previous Krylov basis vectors from Vnext and put coefficients in H and R.
562 
563  // Get a view of all the previous vectors
564  prevind.resize(lclDim);
565  for (int i=0; i<lclDim; i++) { prevind[i] = i; }
566  Vprev = MVT::CloneView(*V_,prevind);
567  Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
568 
569  // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
571  subH = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View,*H_,lclDim,1,0,curDim_ ) );
573  AsubH.append( subH );
574 
575  // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
577  subR = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View,*H_,1,1,lclDim,curDim_ ) );
578 
579  // Project out the previous Krylov vectors and normalize the next vector.
580  rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
581 
582  // Copy over the coefficients to R just in case we run into an error.
583  Teuchos::SerialDenseMatrix<int,ScalarType> subR2( Teuchos::View,R_,lclDim+1,1,0,curDim_ );
584  Teuchos::SerialDenseMatrix<int,ScalarType> subH2( Teuchos::View,*H_,lclDim+1,1,0,curDim_ );
585  subR2.assign(subH2);
586 
587  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1,GCRODRIterOrthoFailure, "Belos::GCRODRIter::iterate(): couldn't generate basis of full rank.");
588 
589  // Update the QR factorization of the upper Hessenberg matrix
590  updateLSQR();
591 
592  // Update basis dim
593  curDim_++;
594  } // end while (statusTest == false)
595  }
596  else { // iterate with recycle space
597  while (stest_->checkStatus(this) != Passed && curDim_+1 <= numBlocks_) {
598  iter_++;
599  int lclDim = curDim_ + 1;
600 
601  // Get the current part of the basis.
602  curind[0] = lclDim;
603  Vnext = MVT::CloneViewNonConst(*V_,curind);
604 
605  // Get a view of the previous vectors.
606  // This is used for orthogonalization and for computing V^H K H
607  curind[0] = curDim_;
608  Vprev = MVT::CloneView(*V_,curind);
609 
610  // Compute the next std::vector in the Krylov basis: Vnext = Op*Vprev
611  lp_->apply(*Vprev,*Vnext);
612  Vprev = Teuchos::null;
613 
614  // First, remove the recycled subspace (C) from Vnext and put coefficients in B.
617  subB = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View,*B_,recycledBlocks_,1,0,curDim_ ) );
619  AsubB.append( subB );
620 
621  // Project out the recycled subspace.
622  ortho_->project( *Vnext, AsubB, C );
623 
624  // Now, remove all previous Krylov basis vectors from Vnext and put coefficients in H and R.
625  // Get a view of all the previous vectors
626  prevind.resize(lclDim);
627  for (int i=0; i<lclDim; i++) { prevind[i] = i; }
628  Vprev = MVT::CloneView(*V_,prevind);
629  Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
630 
631  // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
634  AsubH.append( subH );
635 
636  // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
638 
639  // Project out the previous Krylov vectors and normalize the next vector.
640  rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
641 
642  // Copy over the coefficients to R just in case we run into an error.
643  Teuchos::SerialDenseMatrix<int,ScalarType> subR2( Teuchos::View,R_,lclDim+1,1,0,curDim_ );
644  Teuchos::SerialDenseMatrix<int,ScalarType> subH2( Teuchos::View,*H_,lclDim+1,1,0,curDim_ );
645  subR2.assign(subH2);
646 
647  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1,GCRODRIterOrthoFailure, "Belos::GCRODRIter::iterate(): couldn't generate basis of full rank.");
648 
649  // Update the QR factorization of the upper Hessenberg matrix
650  updateLSQR();
651 
652  // Update basis dim
653  curDim_++;
654 
655  } // end while (statusTest == false)
656  } // end if (U_ == Teuchos::null)
657 
658  }
659 
660 
661  template<class ScalarType, class MV, class OP>
663 
664  int i;
665  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
666 
667  // Get correct dimension based on input "dim"
668  // Remember that ortho failures result in an exit before updateLSQR() is called.
669  // Therefore, it is possible that dim == curDim_.
670  int curDim = curDim_;
671  if ( (dim >= curDim_) && (dim < getMaxSubspaceDim()) )
672  curDim = dim;
673 
675  //
676  // Apply previous transformations and compute new transformation to reduce upper-Hessenberg
677  // system to upper-triangular form.
678  //
679  // QR factorization of Least-Squares system with Givens rotations
680  //
681  for (i=0; i<curDim; i++) {
682  //
683  // Apply previous Givens rotations to new column of Hessenberg matrix
684  //
685  blas.ROT( 1, &R_(i,curDim), 1, &R_(i+1, curDim), 1, &cs_[i], &sn_[i] );
686  }
687  //
688  // Calculate new Givens rotation
689  //
690  blas.ROTG( &R_(curDim,curDim), &R_(curDim+1,curDim), &cs_[curDim], &sn_[curDim] );
691  R_(curDim+1,curDim) = zero;
692  //
693  // Update RHS w/ new transformation
694  //
695  blas.ROT( 1, &z_(curDim), 1, &z_(curDim+1), 1, &cs_[curDim], &sn_[curDim] );
696 
697  } // end updateLSQR()
698 
699 } // end Belos namespace
700 
701 #endif /* BELOS_GCRODR_ITER_HPP */
ScalarType * values() const
Collection of types and exceptions used within the Belos solvers.
void TRSM(ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) 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...
int getBlockSize() const
Get the blocksize to be 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.
Array< T > & append(const T &x)
Class which manages the output and verbosity of the Belos solvers.
int shapeUninitialized(OrdinalType numRows, OrdinalType numCols)
void setSize(int recycledBlocks, int numBlocks)
Set the maximum number of blocks used by the iterative solver and the number of recycled vectors...
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
int getRecycledBlocks() const
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem...
Pure virtual base class for defining the status testing capabilities of Belos.
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed...
int curDim
The current dimension of the reduction.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
GCRODRIterOrthoFailure(const std::string &what_arg)
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
Traits class which defines basic operations on multivectors.
virtual ~GCRODRIter()
Destructor.
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< MV > V
The current Krylov basis.
void setRecycledBlocks(int recycledBlocks)
Set the maximum number of recycled blocks used by the iterative solver.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
GCRODRIterInitFailure(const std::string &what_arg)
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem...
void resetNumIters(int iter=0)
Reset the iteration count.
GCRODRIterInitFailure is thrown when the GCRODRIter object is unable to generate an initial iterate i...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
void initialize()
Initialize the solver with empty data. Calling this method will result in error, as GCRODRIter must b...
Structure to contain pointers to GCRODRIter state variables.
OperatorTraits< ScalarType, MV, OP > OPT
int getNumIters() const
Get the current iteration count.
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
Teuchos::ScalarTraits< ScalarType > SCT
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
GCRODRIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void setBlockSize(int blockSize)
Set the blocksize.
Class which defines basic traits for the operator type.
Teuchos::RCP< MV > C
SCT::magnitudeType MagnitudeType
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:28
Belos&#39;s templated virtual class for providing routines for orthogonalization and orthonormzalition of...
GCRODRIter(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)
GCRODRIter constructor with linear problem, solver utilities, and parameter list of solver options...
Belos header file which uses auto-configuration information to include necessary C++ headers...
SerialDenseMatrix< OrdinalType, ScalarType > & assign(const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
MultiVecTraits< ScalarType, MV > MVT
OrdinalType stride() const
bool isInitialized()
States whether the solver has been initialized or not.
int sizeUninitialized(OrdinalType length_in)
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...

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