Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBlockFGmresIter.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_BLOCK_FGMRES_ITER_HPP
11 #define BELOS_BLOCK_FGMRES_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosGmresIteration.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 
28 #include "Teuchos_BLAS.hpp"
31 #include "Teuchos_ScalarTraits.hpp"
33 #include "Teuchos_TimeMonitor.hpp"
34 
48 namespace Belos {
49 
50 template<class ScalarType, class MV, class OP>
51 class BlockFGmresIter : virtual public GmresIteration<ScalarType,MV,OP> {
52 
53  public:
54 
55  //
56  // Convenience typedefs
57  //
62 
64 
65 
76  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
79  Teuchos::ParameterList &params );
80 
82  virtual ~BlockFGmresIter() {};
84 
85 
87 
88 
110  void iterate();
111 
134 
138  void initialize()
139  {
141  initializeGmres(empty);
142  }
143 
153  state.curDim = curDim_;
154  state.V = V_;
155  state.Z = Z_;
156  state.H = H_;
157  state.R = R_;
158  state.z = z_;
159  return state;
160  }
161 
163 
164 
166 
167 
169  int getNumIters() const { return iter_; }
170 
172  void resetNumIters( int iter = 0 ) { iter_ = iter; }
173 
176  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
177 
179 
185 
187 
190  void updateLSQR( int dim = -1 );
191 
193  int getCurSubspaceDim() const {
194  if (!initialized_) return 0;
195  return curDim_;
196  };
197 
199  int getMaxSubspaceDim() const { return blockSize_*numBlocks_; }
200 
202 
203 
205 
206 
208  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
209 
211  int getBlockSize() const { return blockSize_; }
212 
214  void setBlockSize(int blockSize) { setSize( blockSize, numBlocks_ ); }
215 
217  int getNumBlocks() const { return numBlocks_; }
218 
220  void setNumBlocks(int numBlocks) { setSize( blockSize_, numBlocks ); }
221 
228  void setSize(int blockSize, int numBlocks);
229 
231  bool isInitialized() { return initialized_; }
232 
234 
235  private:
236 
237  //
238  // Internal methods
239  //
241  void setStateSize();
242 
243  //
244  // Classes inputed through constructor that define the linear problem to be solved.
245  //
250 
251  //
252  // Algorithmic parameters
253  //
254  // blockSize_ is the solver block size.
255  // It controls the number of vectors added to the basis on each iteration.
256  int blockSize_;
257  // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
258  int numBlocks_;
259 
260  // Storage for QR factorization of the least squares system.
263 
264  //
265  // Current solver state
266  //
267  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
268  // is capable of running; _initialize is controlled by the initialize() member method
269  // For the implications of the state of initialized_, please see documentation for initialize()
270  bool initialized_;
271 
272  // stateStorageInitialized_ specified that the state storage has be initialized to the current
273  // blockSize_ and numBlocks_. This initialization may be postponed if the linear problem was
274  // generated without the right-hand side or solution vectors.
275  bool stateStorageInitialized_;
276 
277  // Current subspace dimension, and number of iterations performed.
278  int curDim_, iter_;
279 
280  //
281  // State Storage
282  //
283  Teuchos::RCP<MV> V_;
284  Teuchos::RCP<MV> Z_;
285  //
286  // Projected matrices
287  // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
288  //
290  //
291  // QR decomposition of Projected matrices for solving the least squares system HY = B.
292  // R_: Upper triangular reduction of H
293  // z_: Q applied to right-hand side of the least squares system
296 };
297 
299  // Constructor.
300  template<class ScalarType, class MV, class OP>
303  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
306  Teuchos::ParameterList &params ):
307  lp_(problem),
308  om_(printer),
309  stest_(tester),
310  ortho_(ortho),
311  blockSize_(0),
312  numBlocks_(0),
313  initialized_(false),
314  stateStorageInitialized_(false),
315  curDim_(0),
316  iter_(0)
317  {
318  // Get the maximum number of blocks allowed for this Krylov subspace
320  ! params.isParameter ("Num Blocks"), std::invalid_argument,
321  "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
322  const int nb = params.get<int> ("Num Blocks");
323 
324  // Set the block size and allocate data.
325  const int bs = params.get ("Block Size", 1);
326  setSize (bs, nb);
327  }
328 
330  // Set the block size and make necessary adjustments.
331  template <class ScalarType, class MV, class OP>
332  void BlockFGmresIter<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
333  {
334  // This routine only allocates space; it doesn't not perform any computation
335  // any change in size will invalidate the state of the solver.
336 
337  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Belos::BlockFGmresIter::setSize was passed a non-positive argument.");
338  if (blockSize == blockSize_ && numBlocks == numBlocks_) {
339  // do nothing
340  return;
341  }
342 
343  if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
344  stateStorageInitialized_ = false;
345 
346  blockSize_ = blockSize;
347  numBlocks_ = numBlocks;
348 
349  initialized_ = false;
350  curDim_ = 0;
351 
352  // Use the current blockSize_ and numBlocks_ to initialize the state storage.
353  setStateSize();
354 
355  }
356 
358  // Setup the state storage.
359  template <class ScalarType, class MV, class OP>
361  {
362  using Teuchos::RCP;
363  using Teuchos::rcp;
365 
366  if (! stateStorageInitialized_) {
367  // Check if there is any multivector to clone from.
368  RCP<const MV> lhsMV = lp_->getLHS();
369  RCP<const MV> rhsMV = lp_->getRHS();
370  if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
371  stateStorageInitialized_ = false;
372  return;
373  }
374  else {
376  // blockSize*numBlocks dependent
377  //
378  int newsd = blockSize_*(numBlocks_+1);
379 
380  if (blockSize_==1) {
381  cs.resize (newsd);
382  sn.resize (newsd);
383  }
384  else {
385  beta.resize (newsd);
386  }
387 
388  // Initialize the state storage
390  blockSize_ * static_cast<ptrdiff_t> (numBlocks_) > MVT::GetGlobalLength (*rhsMV),
391  std::invalid_argument, "Belos::BlockFGmresIter::setStateSize(): "
392  "Cannot generate a Krylov basis with dimension larger the operator!");
393 
394  // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
395  if (V_ == Teuchos::null) {
396  // Get the multivector that is not null.
397  RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
399  tmp == Teuchos::null, std::invalid_argument,
400  "Belos::BlockFGmresIter::setStateSize(): "
401  "linear problem does not specify multivectors to clone from.");
402  V_ = MVT::Clone (*tmp, newsd);
403  }
404  else {
405  // Generate V_ by cloning itself ONLY if more space is needed.
406  if (MVT::GetNumberVecs (*V_) < newsd) {
407  RCP<const MV> tmp = V_;
408  V_ = MVT::Clone (*tmp, newsd);
409  }
410  }
411 
412  if (Z_ == Teuchos::null) {
413  // Get the multivector that is not null.
414  RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
416  tmp == Teuchos::null, std::invalid_argument,
417  "Belos::BlockFGmresIter::setStateSize(): "
418  "linear problem does not specify multivectors to clone from.");
419  Z_ = MVT::Clone (*tmp, newsd);
420  }
421  else {
422  // Generate Z_ by cloning itself ONLY if more space is needed.
423  if (MVT::GetNumberVecs (*Z_) < newsd) {
424  RCP<const MV> tmp = Z_;
425  Z_ = MVT::Clone (*tmp, newsd);
426  }
427  }
428 
429  // Generate H_ only if it doesn't exist, otherwise resize it.
430  if (H_ == Teuchos::null) {
431  H_ = rcp (new SDM (newsd, newsd-blockSize_));
432  }
433  else {
434  H_->shapeUninitialized (newsd, newsd - blockSize_);
435  }
436 
437  // TODO: Insert logic so that Hessenberg matrix can be saved and reduced matrix is stored in R_
438  //R_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( newsd, newsd-blockSize_ ) );
439  // Generate z_ only if it doesn't exist, otherwise resize it.
440  if (z_ == Teuchos::null) {
441  z_ = rcp (new SDM (newsd, blockSize_));
442  }
443  else {
444  z_->shapeUninitialized (newsd, blockSize_);
445  }
446 
447  // State storage has now been initialized.
448  stateStorageInitialized_ = true;
449  }
450  }
451  }
452 
453 
454  template <class ScalarType, class MV, class OP>
457  {
459 
460  Teuchos::RCP<MV> currentUpdate = Teuchos::null;
461  if (curDim_ == 0) {
462  // If this is the first iteration of the Arnoldi factorization,
463  // then there is no update, so return Teuchos::null.
464  return currentUpdate;
465  }
466  else {
467  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero ();
468  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one ();
470 
471  currentUpdate = MVT::Clone (*Z_, blockSize_);
472 
473  // Make a view and then copy the RHS of the least squares problem. DON'T OVERWRITE IT!
474  SDM y (Teuchos::Copy, *z_, curDim_, blockSize_);
475 
476  // Solve the least squares problem.
478  Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,
479  H_->values (), H_->stride (), y.values (), y.stride ());
480 
481  // Compute the current update.
482  std::vector<int> index (curDim_);
483  for (int i = 0; i < curDim_; ++i) {
484  index[i] = i;
485  }
486  Teuchos::RCP<const MV> Zjp1 = MVT::CloneView (*Z_, index);
487  MVT::MvTimesMatAddMv (one, *Zjp1, y, zero, *currentUpdate);
488  }
489  return currentUpdate;
490  }
491 
492 
493  template <class ScalarType, class MV, class OP>
496  getNativeResiduals (std::vector<MagnitudeType> *norms) const
497  {
498  // NOTE: Make sure the incoming std::vector is the correct size!
499  if (norms != NULL && (int)norms->size() < blockSize_) {
500  norms->resize (blockSize_);
501  }
502 
503  if (norms != NULL) {
505  for (int j = 0; j < blockSize_; ++j) {
506  (*norms)[j] = blas.NRM2 (blockSize_, &(*z_)(curDim_, j), 1);
507  }
508  }
509 
510  // FGmres does not return a residual (multi)vector.
511  return Teuchos::null;
512  }
513 
514 
515  template <class ScalarType, class MV, class OP>
518  {
519  using Teuchos::RCP;
520  using Teuchos::rcp;
521  using std::endl;
524  const ScalarType ZERO = STS::zero ();
525  const ScalarType ONE = STS::one ();
526 
527  // Initialize the state storage if it isn't already.
528  if (! stateStorageInitialized_) {
529  setStateSize ();
530  }
531 
533  ! stateStorageInitialized_, std::invalid_argument,
534  "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
535 
536  // NOTE: In BlockFGmresIter, V and Z are required!!! Inconsistent
537  // multivectors widths and lengths will not be tolerated, and will
538  // be treated with exceptions.
539  const char errstr[] = "Belos::BlockFGmresIter::initialize(): The given "
540  "multivectors must have a consistent length and width.";
541 
542  if (! newstate.V.is_null () && ! newstate.z.is_null ()) {
543 
544  // initialize V_,z_, and curDim_
545 
547  MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_),
548  std::invalid_argument, errstr );
550  MVT::GetNumberVecs(*newstate.V) < blockSize_,
551  std::invalid_argument, errstr );
553  newstate.curDim > blockSize_*(numBlocks_+1),
554  std::invalid_argument, errstr );
555 
556  curDim_ = newstate.curDim;
557  const int lclDim = MVT::GetNumberVecs(*newstate.V);
558 
559  // check size of Z
561  newstate.z->numRows() < curDim_ || newstate.z->numCols() < blockSize_,
562  std::invalid_argument, errstr);
563 
564  // copy basis vectors from newstate into V
565  if (newstate.V != V_) {
566  // only copy over the first block and print a warning.
567  if (curDim_ == 0 && lclDim > blockSize_) {
568  std::ostream& warn = om_->stream (Warnings);
569  warn << "Belos::BlockFGmresIter::initialize(): the solver was "
570  << "initialized with a kernel of " << lclDim << endl
571  << "The block size however is only " << blockSize_ << endl
572  << "The last " << lclDim - blockSize_
573  << " vectors will be discarded." << endl;
574  }
575  std::vector<int> nevind (curDim_ + blockSize_);
576  for (int i = 0; i < curDim_ + blockSize_; ++i) {
577  nevind[i] = i;
578  }
579  RCP<const MV> newV = MVT::CloneView (*newstate.V, nevind);
580  RCP<MV> lclV = MVT::CloneViewNonConst (*V_, nevind);
581  MVT::MvAddMv (ONE, *newV, ZERO, *newV, *lclV);
582 
583  // done with local pointers
584  lclV = Teuchos::null;
585  }
586 
587  // put data into z_, make sure old information is not still hanging around.
588  if (newstate.z != z_) {
589  z_->putScalar();
590  SDM newZ (Teuchos::View, *newstate.z, curDim_ + blockSize_, blockSize_);
591  RCP<SDM> lclz;
592  lclz = rcp (new SDM (Teuchos::View, *z_, curDim_ + blockSize_, blockSize_));
593  lclz->assign (newZ);
594  lclz = Teuchos::null; // done with local pointers
595  }
596  }
597  else {
599  newstate.V == Teuchos::null,std::invalid_argument,
600  "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
601 
603  newstate.z == Teuchos::null,std::invalid_argument,
604  "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
605  }
606 
607  // the solver is initialized
608  initialized_ = true;
609  }
610 
611 
612  template <class ScalarType, class MV, class OP>
614  {
615  using Teuchos::Array;
616  using Teuchos::null;
617  using Teuchos::RCP;
618  using Teuchos::rcp;
619  using Teuchos::View;
621 
622  // Allocate/initialize data structures
623  if (initialized_ == false) {
624  initialize();
625  }
626 
627  // Compute the current search dimension.
628  const int searchDim = blockSize_ * numBlocks_;
629 
630  // Iterate until the status test tells us to stop.
631  // Raise an exception if a computed block is not full rank.
632  while (stest_->checkStatus (this) != Passed && curDim_+blockSize_ <= searchDim) {
633  ++iter_;
634 
635  // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
636  const int lclDim = curDim_ + blockSize_;
637 
638  // Get the current part of the basis.
639  std::vector<int> curind (blockSize_);
640  for (int i = 0; i < blockSize_; ++i) {
641  curind[i] = lclDim + i;
642  }
643  RCP<MV> Vnext = MVT::CloneViewNonConst (*V_, curind);
644 
645  // Get a view of the previous vectors.
646  // This is used for orthogonalization and for computing V^H K H.
647  for (int i = 0; i < blockSize_; ++i) {
648  curind[i] = curDim_ + i;
649  }
650  RCP<const MV> Vprev = MVT::CloneView (*V_, curind);
651  RCP<MV> Znext = MVT::CloneViewNonConst (*Z_, curind);
652 
653  // Compute the next (multi)vector in the Krylov basis: Znext = M*Vprev
654  lp_->applyRightPrec (*Vprev, *Znext);
655  Vprev = null;
656 
657  // Compute the next (multi)vector in the Krylov basis: Vnext = A*Znext
658  lp_->applyOp (*Znext, *Vnext);
659  Znext = null;
660 
661  // Remove all previous Krylov basis vectors from Vnext
662  // Get a view of all the previous vectors
663  std::vector<int> prevind (lclDim);
664  for (int i = 0; i < lclDim; ++i) {
665  prevind[i] = i;
666  }
667  Vprev = MVT::CloneView (*V_, prevind);
668  Array<RCP<const MV> > AVprev (1, Vprev);
669 
670  // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
671  RCP<SDM> subH = rcp (new SDM (View, *H_, lclDim, blockSize_, 0, curDim_));
672  Array<RCP<SDM> > AsubH;
673  AsubH.append (subH);
674 
675  // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
676  RCP<SDM> subR = rcp (new SDM (View, *H_, blockSize_, blockSize_, lclDim, curDim_));
677  const int rank = ortho_->projectAndNormalize (*Vnext, AsubH, subR, AVprev);
679  rank != blockSize_, GmresIterationOrthoFailure,
680  "Belos::BlockFGmresIter::iterate(): After orthogonalization, the new "
681  "basis block does not have full rank. It contains " << blockSize_
682  << " vector" << (blockSize_ != 1 ? "s" : "")
683  << ", but its rank is " << rank << ".");
684 
685  //
686  // V has been extended, and H has been extended.
687  //
688  // Update the QR factorization of the upper Hessenberg matrix
689  //
690  updateLSQR ();
691  //
692  // Update basis dim and release all pointers.
693  //
694  Vnext = null;
695  curDim_ += blockSize_;
696  } // end while (statusTest == false)
697  }
698 
699 
700  template<class ScalarType, class MV, class OP>
702  {
705 
706  const ScalarType zero = STS::zero ();
707  const ScalarType two = (STS::one () + STS::one());
708  ScalarType sigma, mu, vscale, maxelem;
710 
711  // Get correct dimension based on input 'dim'. Remember that
712  // orthogonalization failures result in an exit before
713  // updateLSQR() is called. Therefore, it is possible that dim ==
714  // curDim_.
715  int curDim = curDim_;
716  if (dim >= curDim_ && dim < getMaxSubspaceDim ()) {
717  curDim = dim;
718  }
719 
720  // Apply previous transformations, and compute new transformation
721  // to reduce upper Hessenberg system to upper triangular form.
722  // The type of transformation we use depends the block size. We
723  // use Givens rotations for a block size of 1, and Householder
724  // reflectors otherwise.
725  if (blockSize_ == 1) {
726  // QR factorization of upper Hessenberg matrix using Givens rotations
727  for (int i = 0; i < curDim; ++i) {
728  // Apply previous Givens rotations to new column of Hessenberg matrix
729  blas.ROT (1, &(*H_)(i, curDim), 1, &(*H_)(i+1, curDim), 1, &cs[i], &sn[i]);
730  }
731  // Calculate new Givens rotation
732  blas.ROTG (&(*H_)(curDim, curDim), &(*H_)(curDim+1, curDim), &cs[curDim], &sn[curDim]);
733  (*H_)(curDim+1, curDim) = zero;
734 
735  // Update RHS w/ new transformation
736  blas.ROT (1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim]);
737  }
738  else {
739  // QR factorization of least-squares system using Householder reflectors.
740  for (int j = 0; j < blockSize_; ++j) {
741  // Apply previous Householder reflectors to new block of Hessenberg matrix
742  for (int i = 0; i < curDim + j; ++i) {
743  sigma = blas.DOT (blockSize_, &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
744  sigma += (*H_)(i,curDim+j);
745  sigma *= beta[i];
746  blas.AXPY (blockSize_, ScalarType(-sigma), &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
747  (*H_)(i,curDim+j) -= sigma;
748  }
749 
750  // Compute new Householder reflector
751  const int maxidx = blas.IAMAX (blockSize_+1, &(*H_)(curDim+j,curDim+j), 1);
752  maxelem = (*H_)(curDim + j + maxidx - 1, curDim + j);
753  for (int i = 0; i < blockSize_ + 1; ++i) {
754  (*H_)(curDim+j+i,curDim+j) /= maxelem;
755  }
756  sigma = blas.DOT (blockSize_, &(*H_)(curDim + j + 1, curDim + j), 1,
757  &(*H_)(curDim + j + 1, curDim + j), 1);
758  if (sigma == zero) {
759  beta[curDim + j] = zero;
760  } else {
761  mu = STS::squareroot ((*H_)(curDim+j,curDim+j)*(*H_)(curDim+j,curDim+j)+sigma);
762  if (STS::real ((*H_)(curDim + j, curDim + j)) < STM::zero ()) {
763  vscale = (*H_)(curDim+j,curDim+j) - mu;
764  } else {
765  vscale = -sigma / ((*H_)(curDim+j, curDim+j) + mu);
766  }
767  beta[curDim+j] = two * vscale * vscale / (sigma + vscale*vscale);
768  (*H_)(curDim+j, curDim+j) = maxelem*mu;
769  for (int i = 0; i < blockSize_; ++i) {
770  (*H_)(curDim+j+1+i,curDim+j) /= vscale;
771  }
772  }
773 
774  // Apply new Householder reflector to the right-hand side.
775  for (int i = 0; i < blockSize_; ++i) {
776  sigma = blas.DOT (blockSize_, &(*H_)(curDim+j+1,curDim+j),
777  1, &(*z_)(curDim+j+1,i), 1);
778  sigma += (*z_)(curDim+j,i);
779  sigma *= beta[curDim+j];
780  blas.AXPY (blockSize_, ScalarType(-sigma), &(*H_)(curDim+j+1,curDim+j),
781  1, &(*z_)(curDim+j+1,i), 1);
782  (*z_)(curDim+j,i) -= sigma;
783  }
784  }
785  } // end if (blockSize_ == 1)
786 
787  // If the least-squares problem is updated wrt "dim" then update curDim_.
788  if (dim >= curDim_ && dim < getMaxSubspaceDim ()) {
789  curDim_ = dim + blockSize_;
790  }
791  } // end updateLSQR()
792 
793 } // namespace Belos
794 
795 #endif /* BELOS_BLOCK_FGMRES_ITER_HPP */
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
void iterate()
This method performs block FGmres iterations until the status test indicates the need to stop or an e...
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
void AXPY(const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
void setBlockSize(int blockSize)
Set the blocksize.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
Class which manages the output and verbosity of the Belos solvers.
BlockFGmresIter(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)
BlockFGmresIter constructor with linear problem, solver utilities, and parameter list of solver optio...
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::ScalarTraits< ScalarType > SCT
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
Teuchos::RCP< const MV > V
The current Krylov basis.
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)
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
Declaration of basic traits for the multivector type.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem...
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
OrdinalType IAMAX(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Structure to contain pointers to GmresIteration state variables.
A pure virtual class for defining the status tests for the Belos iterative solvers.
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear prob...
void resetNumIters(int iter=0)
Reset the iteration count.
Class which defines basic traits for the operator type.
int getNumIters() const
Get the current iteration count.
Traits class which defines basic operations on multivectors.
bool isParameter(const std::string &name) const
GmresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const MV > Z
The current preconditioned Krylov basis (only used in flexible GMRES).
int curDim
The current dimension of the reduction.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
This class implements the block flexible GMRES iteration, where a block Krylov subspace is constructe...
ScalarType DOT(const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
virtual ~BlockFGmresIter()
Destructor.
void initializeGmres(GmresIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
bool isInitialized()
States whether the solver has been initialized or not.
MultiVecTraits< ScalarType, MV > MVT
OrdinalType numCols() const
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
Class which defines basic traits for the operator type.
SCT::magnitudeType MagnitudeType
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
Belos&#39;s templated virtual class for providing routines for orthogonalization and orthonormzalition of...
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
OrdinalType numRows() const
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
bool is_null() const

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