Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosLSQRIter.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_LSQR_ITER_HPP
11 #define BELOS_LSQR_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosLSQRIteration.hpp"
20 
21 #include "BelosLinearProblem.hpp"
22 #include "BelosOutputManager.hpp"
23 #include "BelosStatusTest.hpp"
24 #include "BelosOperatorTraits.hpp"
25 #include "BelosMultiVecTraits.hpp"
26 
29 #include "Teuchos_ScalarTraits.hpp"
31 #include "Teuchos_TimeMonitor.hpp"
32 
40 namespace Belos {
41 
42 template<class ScalarType, class MV, class OP>
43 class LSQRIter : virtual public Belos::Iteration<ScalarType,MV,OP> {
44 
45  public:
46 
47  //
48  // Convenience typedefs
49  //
54 
56 
57 
67  Teuchos::ParameterList &params );
68 
69 // If either blocks or reorthogonalization exist, then
70 // const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
71 
72 
74  virtual ~LSQRIter() {};
76 
77 
79 
80 
92  void iterate();
93 
104 
107  void initialize()
108  {
110  initializeLSQR(empty);
111  }
112 
122  state.U = U_; // right Lanczos vector
123  state.V = V_; // left Lanczos vector
124  state.W = W_; // OP * V
125  state.lambda = lambda_;
126  state.resid_norm = resid_norm_;
127  state.frob_mat_norm = frob_mat_norm_;
128  state.mat_cond_num = mat_cond_num_;
129  state.mat_resid_norm = mat_resid_norm_;
130  state.sol_norm = sol_norm_;
131  state.bnorm = bnorm_;
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 { return Teuchos::null; }
150 
152 
154  Teuchos::RCP<MV> getCurrentUpdate() const { return Teuchos::null; }
155 
157 
159 
160 
162  const Belos::LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
163 
165  int getBlockSize() const { return 1; }
166 
167 
169  //This is unique to single vector methods.
170  void setBlockSize(int blockSize) {
171  TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
172  "LSQRIter::setBlockSize(): Cannot use a block size that is not one.");
173  }
174 
176  bool isInitialized() { return initialized_; }
177 
179 
180  private:
181 
182  //
183  // Internal methods
184  //
186  void setStateSize();
187 
188  //
189  // Classes inputed through constructor that define the linear problem to be solved.
190  //
194 // const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
195 
196  //
197  // Current solver state
198  //
199  // initialized_ specifies that the basis vectors have been initialized and the iterate()
200  // routine is capable of running; _initialize is controlled by the initialize() member
201  // method. For the implications of the state of initialized_, please see documentation
202  // for initialize()
203  bool initialized_;
204 
205  // stateStorageInitialized_ specifies that the state storage has been initialized.
206  // This initialization may be postponed if the linear problem was generated without
207  // the right-hand side or solution vectors.
208  bool stateStorageInitialized_;
209 
210  // Current number of iterations performed.
211  int iter_;
212 
213  //
214  // State Storage
215  //
216  //
217  // Bidiagonalization vector
218  Teuchos::RCP<MV> U_;
219  //
220  // Bidiagonalization vector
221  Teuchos::RCP<MV> V_;
222  //
223  // Direction vector
224  Teuchos::RCP<MV> W_;
225  //
226  // Damping value
227  MagnitudeType lambda_;
228  //
229  // Residual norm estimate
230  ScalarType resid_norm_;
231  //
232  // Frobenius norm estimate
233  ScalarType frob_mat_norm_;
234  //
235  // Condition number estimate
236  ScalarType mat_cond_num_;
237  //
238  // A^T*resid norm estimate
239  ScalarType mat_resid_norm_;
240  //
241  // Solution norm estimate
242  ScalarType sol_norm_;
243  //
244  // RHS norm
245  ScalarType bnorm_;
246 
247 };
248 
250  // Constructor.
251 
252 // const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
253 
254  template<class ScalarType, class MV, class OP>
258  Teuchos::ParameterList &params ):
259  lp_(problem),
260  om_(printer),
261  stest_(tester),
262  initialized_(false),
263  stateStorageInitialized_(false),
264  iter_(0),
265  lambda_(params.get<MagnitudeType> ("Lambda")),
266  resid_norm_(0.0),
267  frob_mat_norm_(0.0),
268  mat_cond_num_(0.0),
269  mat_resid_norm_(0.0),
270  sol_norm_(0.0),
271  bnorm_(0.0)
272  {
273  }
274 
276  // Setup the state storage.
277  template <class ScalarType, class MV, class OP>
279  {
280  if (!stateStorageInitialized_) {
281  // Check if there is any multivector to clone from.
282  Teuchos::RCP<const MV> rhsMV = lp_->getInitPrecResVec();
283  Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
284  if (lhsMV == Teuchos::null || rhsMV == Teuchos::null) {
285  stateStorageInitialized_ = false;
286  return;
287  }
288  else {
289  // Initialize the state storage
290  // If the subspace has not been initialized before, generate it
291  // using the LHS and RHS from lp_.
292  if (U_ == Teuchos::null) {
293  // Get the multivectors.
294  TEUCHOS_TEST_FOR_EXCEPTION(rhsMV == Teuchos::null, std::invalid_argument, "LSQRIter::setStateSize(): linear problem does not specify right hand multivector to clone from.");
295  TEUCHOS_TEST_FOR_EXCEPTION(lhsMV == Teuchos::null, std::invalid_argument, "LSQRIter::setStateSize(): linear problem does not specify left hand multivector to clone from.");
296 
297  U_ = MVT::Clone( *rhsMV, 1 ); // LeftPrecond * rhs
298  V_ = MVT::Clone( *lhsMV, 1 ); // zero, overwrittein in
299  W_ = MVT::Clone( *lhsMV, 1 ); // zero, initializeLSQR
300  }
301 
302  // State storage has now been initialized.
303  stateStorageInitialized_ = true;
304  }
305  }
306  }
307 
308 
310  // Initialize this iteration object
311  template <class ScalarType, class MV, class OP>
313  {
314  using Teuchos::RCP;
315 
316  // Initialize the state storage if it isn't already.
317  if (!stateStorageInitialized_)
318  setStateSize();
319 
320  TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
321  "LSQRIter::initialize(): Cannot initialize state storage!");
322 
323  std::string errstr("LSQRIter::initialize(): Specified multivectors must have a consistent length and width.");
324 
325 
326  // Compute initial bidiagonalization vectors and search direction
327  //
328  RCP<const MV> lhsMV = lp_->getLHS(); // contains initial guess,
329 
330  const bool debugSerialLSQR = false;
331 
332  if( debugSerialLSQR )
333  {
334  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> lhsNorm(1);
335  MVT::MvNorm( *lhsMV, lhsNorm );
336  std::cout << "initializeLSQR lhsNorm " << lhsNorm[0] << std::endl;
337  }
338 
339  // LinearProlbem provides right-hand side vectors including RHS CurrRHSVec InitResVec.
340  RCP<const MV> rhsMV = lp_->getInitPrecResVec();
341  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
342  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
343  MVT::MvAddMv( one, *rhsMV, zero, *U_, *U_);
344 
345  RCP<const OP> M_left = lp_->getLeftPrec();
346  RCP<const OP> A = lp_->getOperator();
347  RCP<const OP> M_right = lp_->getRightPrec();
348 
349  if( debugSerialLSQR )
350  {
351  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> rhsNorm(1);
352  MVT::MvNorm( *U_, rhsNorm );
353  std::cout << "initializeLSQR | U_ | : " << rhsNorm[0] << std::endl;
354  //U_->Print(std::cout);
355  }
356 
357  //MVT::MvScale( *V_, zero );
358 
359  // Apply the (conjugate) transpose of the preconditioned operator:
360  //
361  // V := (M_L A M_R)^* U, which means
362  // V := M_R^* (A^* (M_L^* U)).
363  //
364  //OPT::Apply(*(lp_->getOperator()), *U_, *V_, CONJTRANS);
365  if ( M_left.is_null())
366  {
367  OPT::Apply (*A, *U_, *V_, CONJTRANS); // V_ = A' U_
368  //std::cout << "*************** V_ ****************" << std::endl;
369  //V_->Print(std::cout);
370  }
371  else
372  {
373  RCP<MV> tempInRangeOfA = MVT::CloneCopy (*U_);
374  OPT::Apply (*M_left, *U_, *tempInRangeOfA, CONJTRANS);
375  OPT::Apply (*A, *tempInRangeOfA, *V_, CONJTRANS); // V_ = A' LeftPrec' U_
376  //std::cout << "mLeft V_ = " << *V_ << std::endl;
377  }
378  if (! M_right.is_null())
379  {
380  RCP<MV> tempInDomainOfA = MVT::CloneCopy (*V_);
381 
382  OPT::Apply (*M_right, *tempInDomainOfA, *V_, CONJTRANS); // V:= RtPrec' A' LeftPrec' U
383  //std::cout << "mRight V_ = " << *V_ << std::endl;
384  }
385 
386  // W := V (copy the vector)
387  MVT::MvAddMv( one, *V_, zero, *V_, *W_);
388 
389  frob_mat_norm_ = zero; // These are
390  mat_cond_num_ = one; // lower
391  sol_norm_ = zero; // bounds.
392 
393  // The solver is initialized.
394  initialized_ = true;
395  }
396 
397 
399  // Iterate until the status test informs us we should stop.
400  template <class ScalarType, class MV, class OP>
402  {
403  //
404  // Allocate/initialize data structures
405  //
406  if (initialized_ == false) {
407  initialize();
408  }
409 
410  // Create convenience variables for zero and one.
411  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
413  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
414 
415  // Allocate memory for scalars.
416  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> alpha(1);
417  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> beta(1);
418  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> xi(1);
419  // xi is a dumb scalar used for storing inner products.
420  // Eventually SDM will replace the "vectors".
421  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> wnorm2(1);
422  ScalarType rhobar, phibar, cs1, phi, rho, cs, sn, theta, xxnorm = MTzero, common;
423  ScalarType zetabar, sn1, psi, res = zero, ddnorm = zero, gamma, tau;
424  ScalarType anorm2 = zero;
425  ScalarType cs2 = -one, sn2 = zero, gammabar, zeta = zero, delta;
426 
427  // The pair of work vectors AV and AtU are
428  Teuchos::RCP<MV> AV; // used in applying A to V_ and
429  AV = MVT::Clone( *U_, 1);
430  Teuchos::RCP<MV> AtU; // used in applying A^TRANS to U_ respectively.
431  AtU = MVT::Clone( *V_, 1);
432  const bool debugSerialLSQR = false;
433 
434  // Get the current solution vector.
435  Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
436 
437 
438  // Check that the current solution vector only has one column.
439  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*cur_soln_vec) != 1, LSQRIterateFailure,
440  "LSQRIter::iterate(): current linear system has more than one vector!" );
441 
442  // In initializeLSQR among other things V = A' U.
443  // alpha and beta normalize these vectors.
444  MVT::MvNorm( *U_, beta );
445  if( SCT::real(beta[0]) > zero )
446  {
447  MVT::MvScale( *U_, one / beta[0] );
448 
449  //std::cout << "*************** U/beta ****************" << std::endl;
450  //U_->Print(std::cout);
451 
452  MVT::MvScale( *V_, one / beta[0] ); // scale V = A'U to normalize U
453 
454  //std::cout << "*************** V/beta ****************" << std::endl;
455  //V_->Print(std::cout);
456  }
457  MVT::MvNorm( *V_, alpha );
458  if( debugSerialLSQR )
459  {
460  // used to compare with implementations
461  // initializing mat_resid_norm to alpha/beta
462  std::cout << iter_ << " First alpha " << alpha[0] << " beta " << beta[0] << " lambda " << lambda_ << std::endl;
463  }
464  if( SCT::real(alpha[0]) > zero )
465  {
466  MVT::MvScale( *V_, one / alpha[0] ); // V alpha = A' U to normalize V
467  }
468  if( beta[0] * alpha[0] > zero )
469  {
470  MVT::MvScale( *W_, one / (beta[0] * alpha[0]) ); // W = V
471  }
472  else
473  {
474  MVT::MvScale( *W_, zero );
475  }
476 
477  using Teuchos::RCP;
478  RCP<const OP> M_left = lp_->getLeftPrec();
479  RCP<const OP> A = lp_->getOperator();
480  RCP<const OP> M_right = lp_->getRightPrec();
481 
482  rhobar = alpha[0];
483  phibar = beta[0];
484  bnorm_ = beta[0];
485  resid_norm_ = beta[0];
486  mat_resid_norm_ = alpha[0] * beta[0];
487 
488 
490  // Iterate until the status test tells us to stop.
491  //
492  while (stest_->checkStatus(this) != Belos::Passed) {
493  // Increment the iteration
494  iter_++;
495 
496  // Perform the next step of the bidiagonalization.
497  // The next U_ and V_ vectors and scalars alpha and beta satisfy
498  // U_ betaNew := AV - U_ alphaOld ...
499 
500  if ( M_right.is_null() )
501  {
502  OPT::Apply(*A, *V_, *AV); // AV := A * V_
503  }
504  else
505  {
506  RCP<MV> tempInDomainOfA = MVT::CloneCopy (*V_);
507  OPT::Apply (*M_right, *V_, *tempInDomainOfA);
508  OPT::Apply(*A, *tempInDomainOfA, *AV);
509  }
510 
511  if (! M_left.is_null())
512  {
513  RCP<MV> tempInRangeOfA = MVT::CloneCopy (*AV);
514  OPT::Apply (*M_left, *tempInRangeOfA, *AV); // AV may change
515  }
516 
517 
518  if ( !( M_left.is_null() && M_right.is_null() )
519  && debugSerialLSQR && iter_ == 1)
520  {
521  // In practice, LSQR may reveal bugs in transposed preconditioners.
522  // This is the test that catches this type of bug.
523  // 1. confirm that V alpha = A' U
524 
525  if (! M_left.is_null())
526  {
527  RCP<MV> tempInRangeOfA = MVT::CloneCopy (*U_);
528  OPT::Apply (*M_left, *U_, *tempInRangeOfA, CONJTRANS);
529  OPT::Apply (*A, *tempInRangeOfA, *AtU, CONJTRANS); // AtU = B'L'U
530  }
531  else
532  {
533  OPT::Apply (*A, *U_, *AtU, CONJTRANS); // AtU = B'U
534  }
535  if ( !( M_right.is_null() ) )
536  {
537  RCP<MV> tempInDomainOfA = MVT::CloneCopy (*AtU);
538  OPT::Apply (*M_right, *tempInDomainOfA, *AtU, CONJTRANS); // AtU := R' AtU
539  }
540 
541  MVT::MvAddMv( one, *AtU, -alpha[0], *V_, *AtU );
542  MVT::MvNorm( *AtU, xi );
543  std::cout << "| V alpha - A' u |= " << xi[0] << std::endl;
544  // 2. confirm that U is a unit vector
546  MVT::MvTransMv( one, *U_, *U_, uotuo );
547  std::cout << "<U, U> = " << printMat(uotuo) << std::endl;
548  // 3. print alpha = <V, A'U>
549  std::cout << "alpha = " << alpha[0] << std::endl;
550  // 4. compute < AV, U> which ought to be alpha
552  MVT::MvTransMv( one, *AV, *U_, utav );
553  std::cout << "<AV, U> = alpha = " << printMat(utav) << std::endl;
554  }
555 
556  MVT::MvAddMv( one, *AV, -alpha[0], *U_, *U_ ); // uNew := Av - uOld alphaOld
557  MVT::MvNorm( *U_, beta);
558 
559  anorm2 += alpha[0]*alpha[0] + beta[0]*beta[0] + lambda_*lambda_;
560 
561 
562  if ( SCT::real(beta[0]) > zero )
563  {
564 
565  MVT::MvScale( *U_, one / beta[0] );
566 
567  if (M_left.is_null())
568  { // ... and V_ alphaNew := AtU - V_ betaNew
569  OPT::Apply(*A, *U_, *AtU, CONJTRANS);
570  }
571  else
572  {
573  RCP<MV> tempInRangeOfA = MVT::CloneCopy (*U_);
574  OPT::Apply (*M_left, *U_, *tempInRangeOfA, CONJTRANS);
575  OPT::Apply(*A, *tempInRangeOfA, *AtU, CONJTRANS);
576  }
577  if (! M_right.is_null())
578  {
579  RCP<MV> tempInDomainOfA = MVT::CloneCopy (*AtU);
580  OPT::Apply (*M_right, *tempInDomainOfA, *AtU, CONJTRANS); // AtU may change
581  }
582 
583  MVT::MvAddMv( one, *AtU, -beta[0], *V_, *V_ );
584  MVT::MvNorm( *V_, alpha );
585  }
586  else // beta = 0
587  {
588  alpha[0] = zero;
589  }
590 
591  if ( SCT::real(alpha[0]) > zero )
592  {
593  MVT::MvScale( *V_, one / alpha[0] );
594  }
595 
596  // Use a plane rotation to eliminate the damping parameter.
597  // This alters the diagonal (rhobar) of the lower-bidiagonal matrix.
598  common = Teuchos::ScalarTraits< ScalarType >::squareroot(rhobar*rhobar + lambda_*lambda_);
599  cs1 = rhobar / common;
600  sn1 = lambda_ / common;
601  psi = sn1 * phibar;
602  phibar = cs1 * phibar;
603 
604  // Use a plane rotation to eliminate the subdiagonal element (beta)
605  // of the lower-bidiagonal matrix, giving an upper-bidiagonal matrix.
606  rho = Teuchos::ScalarTraits< ScalarType >::squareroot(rhobar*rhobar + lambda_*lambda_ + beta[0]*beta[0]);
607  cs = common / rho;
608  sn = beta[0] / rho;
609  theta = sn * alpha[0];
610  rhobar = -cs * alpha[0];
611  phi = cs * phibar;
612  phibar = sn * phibar; // If beta vanishes, so do sn, theta, phibar and eventually resid_norm
613  tau = sn * phi;
614 
615  delta = sn2 * rho;
616  gammabar = -cs2 * rho;
617  zetabar = (phi - delta*zeta) / gammabar;
618  sol_norm_ = Teuchos::ScalarTraits< ScalarType >::squareroot(xxnorm + zetabar*zetabar);
619  gamma = Teuchos::ScalarTraits< ScalarType >::squareroot(gammabar*gammabar + theta*theta);
620  cs2 = gammabar / gamma;
621  sn2 = theta / gamma;
622  zeta = (phi - delta*zeta) / gamma;
623  xxnorm += zeta*zeta;
624 
625  //The next task may be addressed by some form of lp_->updateSolution.
626  if ( M_right.is_null())
627  {
628  MVT::MvAddMv( phi / rho, *W_, one, *cur_soln_vec, *cur_soln_vec);
629  }
630  else
631  {
632  RCP<MV> tempInDomainOfA = MVT::CloneCopy (*W_);
633  OPT::Apply (*M_right, *W_, *tempInDomainOfA);
634  MVT::MvAddMv( phi / rho, *tempInDomainOfA, one, *cur_soln_vec, *cur_soln_vec);
635  }
636 
637  MVT::MvNorm( *W_, wnorm2 );
638  ddnorm += (one / rho)*(one / rho) * wnorm2[0]*wnorm2[0];
639  MVT::MvAddMv( one, *V_, -theta / rho, *W_, *W_ );
640 
641  frob_mat_norm_ = Teuchos::ScalarTraits< ScalarType >::squareroot(anorm2);
642  mat_cond_num_ = frob_mat_norm_ * Teuchos::ScalarTraits< ScalarType >::squareroot(ddnorm);
643  res+= psi*psi;
644  resid_norm_ = Teuchos::ScalarTraits< ScalarType >::squareroot(phibar*phibar + res);
645  mat_resid_norm_ = alpha[0] * Teuchos::ScalarTraits< ScalarType >::magnitude(tau);
646 
647  } // end while (sTest_->checkStatus(this) != Passed)
648  } // iterate()
649 
650 } // end Belos namespace
651 
652 #endif /* BELOS_LSQR_ITER_HPP */
Teuchos::RCP< const MV > U
Bidiagonalization vector.
Collection of types and exceptions used within the Belos solvers.
IterationState contains the data that defines the state of the LSQR solver at any given time...
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.
ScalarType resid_norm
The current residual norm.
Teuchos::ScalarTraits< ScalarType >::magnitudeType lambda
The damping value.
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
Get the norms of the residuals native to the solver.
SerialBandDenseMatrixPrinter< OrdinalType, ScalarType > printMat(const SerialBandDenseMatrix< OrdinalType, ScalarType > &obj)
Pure virtual base class for defining the status testing capabilities of Belos.
static T squareroot(T x)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
void initializeLSQR(LSQRIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, completing the initial state.
ScalarType sol_norm
An estimate of the norm of the solution.
ScalarType mat_cond_num
An approximation to the condition number of A.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Belos::OperatorTraits< ScalarType, MV, OP > OPT
Class which defines basic traits for the operator type.
Traits class which defines basic operations on multivectors.
void iterate()
This method performs LSQR iterations until the status test indicates the need to stop or an error occ...
ScalarType frob_mat_norm
An approximation to the Frobenius norm of A.
ScalarType mat_resid_norm
An estimate of the norm of A^T*resid.
Teuchos::RCP< const MV > W
The search direction vector.
Teuchos::RCP< const MV > V
Bidiagonalization vector.
virtual ~LSQRIter()
Destructor.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getNumIters() const
Get the current iteration count.
LSQRIter(const Teuchos::RCP< Belos::LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Belos::OutputManager< ScalarType > > &printer, const Teuchos::RCP< Belos::StatusTest< ScalarType, MV, OP > > &tester, Teuchos::ParameterList &params)
LSQRIter constructor with linear problem, solver utilities, and parameter list of solver options...
Teuchos::ScalarTraits< ScalarType > SCT
static magnitudeType magnitude(T a)
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Belos::MultiVecTraits< ScalarType, MV > MVT
SCT::magnitudeType MagnitudeType
LSQRIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
void initialize()
The solver is initialized using initializeLSQR.
Class which defines basic traits for the operator type.
Implementation of the LSQR iteration.
ScalarType bnorm
The norm of the RHS vector b.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
bool isInitialized()
States whether the solver has been initialized or not.
Structure to contain pointers to LSQRIteration state variables, ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
const Belos::LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
LSQRIterateFailure is thrown when the LSQRIteration object is unable to compute the next iterate in t...
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver to solve this linear problem.

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