Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosGmresPolyOp.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 
42 #ifndef BELOS_GMRESPOLYOP_HPP
43 #define BELOS_GMRESPOLYOP_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosOperator.hpp"
53 #include "BelosMultiVec.hpp"
54 #include "BelosOperatorTraits.hpp"
55 #include "BelosMultiVecTraits.hpp"
56 #include "BelosLinearProblem.hpp"
57 
58 #include "BelosGmresIteration.hpp"
59 #include "BelosBlockGmresIter.hpp"
60 
64 
68 #include "BelosStatusTestCombo.hpp"
70 
71 #include "BelosOutputManager.hpp"
72 
73 #include "Teuchos_BLAS.hpp"
74 #include "Teuchos_LAPACK.hpp"
75 #include "Teuchos_as.hpp"
76 #include "Teuchos_RCP.hpp"
81 
82 #ifdef BELOS_TEUCHOS_TIME_MONITOR
83  #include "Teuchos_TimeMonitor.hpp"
84 #endif // BELOS_TEUCHOS_TIME_MONITOR
85 
86 namespace Belos {
87 
94  class GmresPolyOpOrthoFailure : public BelosError {public:
95  GmresPolyOpOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
96  {}};
97 
98  // Create a shell class for the MV, inherited off MultiVec<> that will operate with the GmresPolyOp.
99  template <class ScalarType, class MV>
100  class GmresPolyMv : public MultiVec< ScalarType >
101  {
102  public:
103 
104  GmresPolyMv ( const Teuchos::RCP<MV>& mv_in )
105  : mv_(mv_in)
106  {}
108  {
109  mv_ = Teuchos::rcp_const_cast<MV>( mv_in );
110  }
111  Teuchos::RCP<MV> getMV() { return mv_; }
112  Teuchos::RCP<const MV> getConstMV() const { return mv_; }
113 
114  GmresPolyMv * Clone ( const int numvecs ) const
115  {
116  GmresPolyMv * newMV = new GmresPolyMv( MVT::Clone( *mv_, numvecs ) );
117  return newMV;
118  }
119  GmresPolyMv * CloneCopy () const
120  {
121  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneCopy( *mv_ ) );
122  return newMV;
123  }
124  GmresPolyMv * CloneCopy ( const std::vector<int>& index ) const
125  {
126  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneCopy( *mv_, index ) );
127  return newMV;
128  }
129  GmresPolyMv * CloneViewNonConst ( const std::vector<int>& index )
130  {
131  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneViewNonConst( *mv_, index ) );
132  return newMV;
133  }
134  const GmresPolyMv * CloneView ( const std::vector<int>& index ) const
135  {
136  const GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneView( *mv_, index ) );
137  return newMV;
138  }
139  ptrdiff_t GetGlobalLength () const { return MVT::GetGlobalLength( *mv_ ); }
140  int GetNumberVecs () const { return MVT::GetNumberVecs( *mv_ ); }
141  void MvTimesMatAddMv (const ScalarType alpha,
142  const MultiVec<ScalarType>& A,
143  const Teuchos::SerialDenseMatrix<int,ScalarType>& B, const ScalarType beta)
144  {
145  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
146  MVT::MvTimesMatAddMv( alpha, *(A_in.getConstMV()), B, beta, *mv_ );
147  }
148  void MvAddMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, const ScalarType beta, const MultiVec<ScalarType>& B )
149  {
150  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
151  const GmresPolyMv<ScalarType,MV>& B_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(B);
152  MVT::MvAddMv( alpha, *(A_in.getConstMV()), beta, *(B_in.getConstMV()), *mv_ );
153  }
154  void MvScale ( const ScalarType alpha ) { MVT::MvScale( *mv_, alpha ); }
155  void MvScale ( const std::vector<ScalarType>& alpha ) { MVT::MvScale( *mv_, alpha ); }
156  void MvTransMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, Teuchos::SerialDenseMatrix<int,ScalarType>& B) const
157  {
158  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
159  MVT::MvTransMv( alpha, *(A_in.getConstMV()), *mv_, B );
160  }
161  void MvDot ( const MultiVec<ScalarType>& A, std::vector<ScalarType>& b ) const
162  {
163  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
164  MVT::MvDot( *(A_in.getConstMV()), *mv_, b );
165  }
166  void MvNorm ( std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec, NormType type = TwoNorm ) const
167  {
168  MVT::MvNorm( *mv_, normvec, type );
169  }
170  void SetBlock ( const MultiVec<ScalarType>& A, const std::vector<int>& index )
171  {
172  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
173  MVT::SetBlock( *(A_in.getConstMV()), index, *mv_ );
174  }
175  void MvRandom () { MVT::MvRandom( *mv_ ); }
176  void MvInit ( const ScalarType alpha ) { MVT::MvInit( *mv_, alpha ); }
177  void MvPrint ( std::ostream& os ) const { MVT::MvPrint( *mv_, os ); }
178 
179  private:
180 
181  typedef MultiVecTraits<ScalarType,MV> MVT;
182 
183  Teuchos::RCP<MV> mv_;
184 
185  };
186 
197  template <class ScalarType, class MV, class OP>
198  class GmresPolyOp : public Operator<ScalarType> {
199  public:
200 
202 
203 
206  const Teuchos::RCP<Teuchos::ParameterList>& params_in
207  )
208  : problem_(problem_in),
209  params_(params_in),
210  LP_(problem_in->getLeftPrec()),
211  RP_(problem_in->getRightPrec())
212  {
213  setParameters( params_ );
214 
215  polyUpdateLabel_ = label_ + ": Hybrid Gmres: Vector Update";
216 #ifdef BELOS_TEUCHOS_TIME_MONITOR
217  timerPolyUpdate_ = Teuchos::TimeMonitor::getNewCounter(polyUpdateLabel_);
218 #endif // BELOS_TEUCHOS_TIME_MONITOR
219 
220  if (polyType_ == "Arnoldi" || polyType_=="Roots")
222  else if (polyType_ == "Gmres")
224  else
225  TEUCHOS_TEST_FOR_EXCEPTION(polyType_!="Arnoldi"&&polyType_!="Gmres"&&polyType_!="Roots",std::invalid_argument,
226  "Belos::GmresPolyOp: \"Polynomial Type\" must be either \"Arnoldi\", \"Gmres\", or \"Roots\".");
227  }
228 
231  : problem_(problem_in)
232  {
233  // If dimension is zero, it will just apply the operator from problem_in in the Apply method.
234  dim_ = 0;
235  }
236 
238  virtual ~GmresPolyOp() {};
240 
242 
243 
245  void setParameters( const Teuchos::RCP<Teuchos::ParameterList>& params_in );
247 
249 
250 
254  void generateArnoldiPoly();
255 
259  void generateGmresPoly();
260 
262 
264 
265 
271  void ApplyPoly ( const MV& x, MV& y ) const;
272  void ApplyArnoldiPoly ( const MV& x, MV& y ) const;
273  void ApplyGmresPoly ( const MV& x, MV& y ) const;
274  void ApplyRootsPoly ( const MV& x, MV& y ) const;
275 
279  void Apply ( const MultiVec<ScalarType>& x, MultiVec<ScalarType>& y, ETrans /* trans */=NOTRANS ) const
280  {
281  const GmresPolyMv<ScalarType,MV>& x_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(x);
282  GmresPolyMv<ScalarType,MV>& y_in = dynamic_cast<GmresPolyMv<ScalarType,MV>&>(y);
283  ApplyPoly( *(x_in.getConstMV()), *(y_in.getMV()) );
284  }
285 
286  int polyDegree() const { return dim_; }
287 
288  private:
289 
290 #ifdef BELOS_TEUCHOS_TIME_MONITOR
291  Teuchos::RCP<Teuchos::Time> timerPolyUpdate_;
292 #endif // BELOS_TEUCHOS_TIME_MONITOR
293  std::string polyUpdateLabel_;
294 
295  typedef int OT; //Ordinal type
296  typedef MultiVecTraits<ScalarType,MV> MVT;
298  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
300 
301  // Default polynomial parameters
302  static constexpr int maxDegree_default_ = 25;
303  static constexpr int verbosity_default_ = Belos::Errors;
304  static constexpr bool randomRHS_default_ = true;
305  static constexpr const char * label_default_ = "Belos";
306  static constexpr const char * polyType_default_ = "Roots";
307  static constexpr const char * orthoType_default_ = "DGKS";
308  static constexpr std::ostream * outputStream_default_ = &std::cout;
309  static constexpr bool damp_default_ = false;
310  static constexpr bool addRoots_default_ = true;
311 
312  // Variables for generating the polynomial
315  Teuchos::RCP<const OP> LP_, RP_;
316 
317  // Output manager.
319  Teuchos::RCP<std::ostream> outputStream_ = Teuchos::rcp(outputStream_default_,false);
320 
321  // Orthogonalization manager.
323 
324  // Current polynomial parameters
325  MagnitudeType polyTol_ = DefaultSolverParameters::polyTol;
326  int maxDegree_ = maxDegree_default_;
327  int verbosity_ = verbosity_default_;
328  bool randomRHS_ = randomRHS_default_;
329  std::string label_ = label_default_;
330  std::string polyType_ = polyType_default_;
331  std::string orthoType_ = orthoType_default_;
332  int dim_ = 0;
333  bool damp_ = damp_default_;
334  bool addRoots_ = addRoots_default_;
335 
336  // Variables for Arnoldi polynomial
337  mutable Teuchos::RCP<MV> V_, wL_, wR_;
340 
341  // Variables for Gmres polynomial;
342  bool autoDeg = false;
344 
345  // Variables for Roots polynomial:
347 
348  // Modified Leja sorting function. Takes a serial dense matrix of M harmonic Ritz values and an index
349  // of values from 0 to M. Returns the sorted values and sorted index, similar to Matlab.
350  void SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index) const ;
351 
352  //Function determines whether added roots are needed and adds them if option is turned on.
353  void ComputeAddedRoots();
354  };
355 
356  template <class ScalarType, class MV, class OP>
358  {
359  // Check which Gmres polynomial to use
360  if (params_in->isParameter("Polynomial Type")) {
361  polyType_ = params_in->get("Polynomial Type", polyType_default_);
362  }
363 
364  // Check for polynomial convergence tolerance
365  if (params_in->isParameter("Polynomial Tolerance")) {
366  if (params_in->isType<MagnitudeType> ("Polynomial Tolerance")) {
367  polyTol_ = params_in->get ("Polynomial Tolerance",
368  static_cast<MagnitudeType> (DefaultSolverParameters::polyTol));
369  }
370  else {
371  polyTol_ = params_in->get ("Polynomial Tolerance", DefaultSolverParameters::polyTol);
372  }
373  }
374 
375  // Check for maximum polynomial degree
376  if (params_in->isParameter("Maximum Degree")) {
377  maxDegree_ = params_in->get("Maximum Degree", maxDegree_default_);
378  }
379 
380  // Check for maximum polynomial degree
381  if (params_in->isParameter("Random RHS")) {
382  randomRHS_ = params_in->get("Random RHS", randomRHS_default_);
383  }
384 
385  // Check for a change in verbosity level
386  if (params_in->isParameter("Verbosity")) {
387  if (Teuchos::isParameterType<int>(*params_in,"Verbosity")) {
388  verbosity_ = params_in->get("Verbosity", verbosity_default_);
389  }
390  else {
391  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params_in,"Verbosity");
392  }
393  }
394 
395  if (params_in->isParameter("Orthogonalization")) {
396  orthoType_ = params_in->get("Orthogonalization",orthoType_default_);
397  TEUCHOS_TEST_FOR_EXCEPTION( orthoType_ != "DGKS" && orthoType_ != "ICGS" && orthoType_ != "IMGS",
398  std::invalid_argument,
399  "Belos::GmresPolyOp: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
400  }
401 
402  // Check for timer label
403  if (params_in->isParameter("Timer Label")) {
404  label_ = params_in->get("Timer Label", label_default_);
405  }
406 
407  // Output stream
408  if (params_in->isParameter("Output Stream")) {
409  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params_in,"Output Stream");
410  }
411 
412  // Check for damped polynomial
413  if (params_in->isParameter("Damped Poly")) {
414  damp_ = params_in->get("Damped Poly", damp_default_);
415  }
416 
417  // Check for root-adding
418  if (params_in->isParameter("Add Roots")) {
419  addRoots_ = params_in->get("Add Roots", addRoots_default_);
420  }
421  }
422 
423  template <class ScalarType, class MV, class OP>
425  {
426  Teuchos::RCP< MV > V = MVT::Clone( *problem_->getRHS(), maxDegree_+2 );
427 
428  //Make power basis:
429  std::vector<int> index(1,0);
430  Teuchos::RCP< MV > V0 = MVT::CloneViewNonConst(*V, index);
431  if (randomRHS_)
432  MVT::MvRandom( *V0 );
433  else
434  MVT::Assign( *problem_->getRHS(), *V0 );
435 
436  if ( !LP_.is_null() ) {
437  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
438  problem_->applyLeftPrec( *Vtemp, *V0);
439  }
440  if ( damp_ ) {
441  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
442  problem_->apply( *Vtemp, *V0);
443  }
444 
445  for(int i=0; i< maxDegree_+1; i++)
446  {
447  index[0] = i;
448  Teuchos::RCP< const MV > Vi = MVT::CloneView(*V, index);
449  index[0] = i+1;
450  Teuchos::RCP< MV > Vip1 = MVT::CloneViewNonConst(*V, index);
451  problem_->apply( *Vi, *Vip1);
452  }
453 
454  //Consider AV:
455  Teuchos::Range1D range( 1, maxDegree_+1);
456  Teuchos::RCP< const MV > AV = MVT::CloneView( *V, range);
457 
458  //Make lhs (AV)^T(AV)
459  Teuchos::SerialDenseMatrix< OT, ScalarType > AVtransAV( maxDegree_+1, maxDegree_+1);
460  MVT::MvTransMv( SCT::one(), *AV, *AV, AVtransAV);
461  //This process adds pDeg*pDeg + pDeg inner products that aren't in the final count.
462 
464  int infoInt;
465  bool status = true; //Keep adjusting poly deg when true.
466 
467  dim_ = maxDegree_;
469  while( status && dim_ >= 1)
470  {
471  Teuchos::SerialDenseMatrix< OT, ScalarType > lhstemp(Teuchos::Copy, AVtransAV, dim_+1, dim_+1);
472  lapack.POTRF( 'U', dim_+1, lhstemp.values(), lhstemp.stride(), &infoInt);
473 
474  if(autoDeg == false)
475  {
476  status = false;
477  if(infoInt != 0)
478  {
479  std::cout << "BelosGmresPolyOp.hpp: LAPACK POTRF was not successful!!" << std::endl;
480  std::cout << "Error code: " << infoInt << std::endl;
481  }
482  }
483  else
484  {
485  if(infoInt != 0)
486  {//Had bad factor. Reduce poly degree.
487  dim_--;
488  }
489  else
490  {
491  status = false;
492  }
493  }
494  if(status == false)
495  {
496  lhs = lhstemp;
497  }
498  }
499  if(dim_ == 0)
500  {
501  pCoeff_.shape( 1, 1);
502  pCoeff_(0,0) = SCT::one();
503  std::cout << "Poly Degree is zero. No preconditioner created." << std::endl;
504  }
505  else
506  {
507  pCoeff_.shape( dim_+1, 1);
508  //Get correct submatrix of AV:
509  Teuchos::Range1D rangeSub( 1, dim_+1);
510  Teuchos::RCP< const MV > AVsub = MVT::CloneView( *V, rangeSub);
511 
512  //Compute rhs (AV)^T V0
513  MVT::MvTransMv( SCT::one(), *AVsub, *V0, pCoeff_);
514  lapack.POTRS( 'U', dim_+1, 1, lhs.values(), lhs.stride(), pCoeff_.values(), pCoeff_.stride(), &infoInt);
515  if(infoInt != 0)
516  {
517  std::cout << "BelosGmresPolyOp.hpp: LAPACK POTRS was not successful!!" << std::endl;
518  std::cout << "Error code: " << infoInt << std::endl;
519  }
520  }
521  }
522 
523  template <class ScalarType, class MV, class OP>
525  {
526  std::string polyLabel = label_ + ": GmresPolyOp creation";
527 
528  // Create a copy of the linear problem that has a zero initial guess and random RHS.
529  std::vector<int> idx(1,0);
530  Teuchos::RCP<MV> newX = MVT::Clone( *(problem_->getLHS()), 1 );
531  Teuchos::RCP<MV> newB = MVT::Clone( *(problem_->getRHS()), 1 );
532  MVT::MvInit( *newX, SCT::zero() );
533  if (randomRHS_) {
534  MVT::MvRandom( *newB );
535  }
536  else {
537  MVT::Assign( *(MVT::CloneView(*(problem_->getRHS()), idx)), *newB );
538  }
540  Teuchos::rcp( new LinearProblem<ScalarType,MV,OP>( problem_->getOperator(), newX, newB ) );
541  newProblem->setInitResVec( newB );
542  newProblem->setLeftPrec( problem_->getLeftPrec() );
543  newProblem->setRightPrec( problem_->getRightPrec() );
544  newProblem->setLabel(polyLabel);
545  newProblem->setProblem();
546  newProblem->setLSIndex( idx );
547 
548  // Create a parameter list for the GMRES iteration.
549  Teuchos::ParameterList polyList;
550 
551  // Tell the block solver that the block size is one.
552  polyList.set("Num Blocks",maxDegree_);
553  polyList.set("Block Size",1);
554  polyList.set("Keep Hessenberg", true);
555 
556  // Create orthogonalization manager if we need to.
557  if (ortho_.is_null()) {
558  params_->set("Orthogonalization", orthoType_);
559  if (orthoType_=="DGKS") {
560  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( polyLabel ) );
561  }
562  else if (orthoType_=="ICGS") {
563  ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( polyLabel ) );
564  }
565  else if (orthoType_=="IMGS") {
566  ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( polyLabel ) );
567  }
568  else {
569  TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::invalid_argument,
570  "Belos::GmresPolyOp(): Invalid orthogonalization type.");
571  }
572  }
573 
574  // Create output manager.
575  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
576 
577  // Create a simple status test that either reaches the relative residual tolerance or maximum polynomial size.
580 
581  // Implicit residual test, using the native residual to determine if convergence was achieved.
584  convTst->defineScaleForm( convertStringToScaleType("Norm of RHS"), Belos::TwoNorm );
585 
586  // Convergence test that stops the iteration when either are satisfied.
589 
590  // Create Gmres iteration object to perform one cycle of Gmres.
592  gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(newProblem,printer_,polyTest,ortho_,polyList) );
593 
594  // Create the first block in the current Krylov basis (residual).
595  Teuchos::RCP<MV> V_0 = MVT::CloneCopy( *newB );
596  if ( !LP_.is_null() )
597  newProblem->applyLeftPrec( *newB, *V_0 );
598  if ( damp_ )
599  {
600  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V_0);
601  newProblem->apply( *Vtemp, *V_0 );
602  }
603 
604  // Get a matrix to hold the orthonormalization coefficients.
605  r0_.resize(1);
606 
607  // Orthonormalize the new V_0
608  int rank = ortho_->normalize( *V_0, Teuchos::rcpFromRef(r0_) );
610  "Belos::GmresPolyOp::generateArnoldiPoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
611 
612  // Set the new state and initialize the solver.
614  newstate.V = V_0;
615  newstate.z = Teuchos::rcpFromRef( r0_);
616  newstate.curDim = 0;
617  gmres_iter->initializeGmres(newstate);
618 
619  // Perform Gmres iteration
620  try {
621  gmres_iter->iterate();
622  }
623  catch (GmresIterationOrthoFailure e) {
624  // Try to recover the most recent least-squares solution
625  gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
626  }
627  catch (std::exception e) {
628  using std::endl;
629  printer_->stream(Errors) << "Error! Caught exception in BlockGmresIter::iterate() at iteration "
630  << gmres_iter->getNumIters() << endl << e.what () << endl;
631  throw;
632  }
633 
634  // Get the solution for this polynomial, use in comparison below
635  Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
636 
637  // Record polynomial info, get current GMRES state
638  GmresIterationState<ScalarType,MV> gmresState = gmres_iter->getState();
639 
640  // If the polynomial has no dimension, the tolerance is too low, return false
641  dim_ = gmresState.curDim;
642  if (dim_ == 0) {
643  return;
644  }
645  if(polyType_ == "Arnoldi"){
646  // Make a view and then copy the RHS of the least squares problem.
647  //
648  y_ = Teuchos::SerialDenseMatrix<OT,ScalarType>( Teuchos::Copy, *gmresState.z, dim_, 1 );
649  H_ = *gmresState.H;
650 
651  //
652  // Solve the least squares problem.
653  //
656  Teuchos::NON_UNIT_DIAG, dim_, 1, SCT::one(),
657  gmresState.R->values(), gmresState.R->stride(),
658  y_.values(), y_.stride() );
659  }
660  else{ //Generate Roots Poly
661  //Find Harmonic Ritz Values to use as polynomial roots:
662 
663  //Copy of square H used to find poly roots:
664  H_ = Teuchos::SerialDenseMatrix<OT,ScalarType>(Teuchos::Copy, *gmresState.H, dim_, dim_);
665  //Zero out below subdiagonal of H:
666  for(int i=0; i <= dim_-3; i++) {
667  for(int k=i+2; k <= dim_-1; k++) {
668  H_(k,i) = SCT::zero();
669  }
670  }
671  //Extra copy of H because equilibrate changes the matrix:
672  Teuchos::SerialDenseMatrix<OT,ScalarType> Htemp (Teuchos::Copy, H_, dim_, dim_);
673 
674  //View the m+1,m element and last col of H:
675  ScalarType Hlast = (*gmresState.H)(dim_,dim_-1);
676  Teuchos::SerialDenseMatrix<OT,ScalarType> HlastCol (Teuchos::View, H_, dim_, 1, 0, dim_-1);
677 
678  //Set up linear system for H^{-*}e_m:
680  E.putScalar(SCT::zero());
681  E(dim_-1,0) = SCT::one();
682 
684  HSolver.setMatrix( Teuchos::rcpFromRef(Htemp));
686  HSolver.setVectors( Teuchos::rcpFromRef(F), Teuchos::rcpFromRef(E));
687  HSolver.factorWithEquilibration( true );
688 
689  //Factor matrix and solve for F = H^{-*}e_m:
690  int info = 0;
691  info = HSolver.factor();
692  if(info != 0){
693  std::cout << "Hsolver factor: info = " << info << std::endl;
694  }
695  info = HSolver.solve();
696  if(info != 0){
697  std::cout << "Hsolver solve : info = " << info << std::endl;
698  }
699 
700  //Scale F and adjust H for Harmonic Ritz value eigenproblem:
701  F.scale(Hlast*Hlast);
702  HlastCol += F;
703 
704  //Set up for eigenvalue problem to get Harmonic Ritz Values:
706  theta_.shape(dim_,2);//1st col for real part, 2nd col for imaginary
707 
708  const int ldv = 1;
709  ScalarType* vlr = 0;
710 
711  // Size of workspace and workspace for DGEEV
712  int lwork = -1;
713  std::vector<ScalarType> work(1);
714  std::vector<MagnitudeType> rwork(2*dim_);
715 
716  //Find workspace size for DGEEV:
717  lapack.GEEV('N','N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
718  lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
719  work.resize( lwork );
720  // Solve for Harmonic Ritz Values:
721  lapack.GEEV('N','N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
722 
723  if(info != 0){
724  std::cout << "GEEV solve : info = " << info << std::endl;
725  }
726 
727  // Set index for sort function and sort Harmonic Ritz Values:
728  std::vector<int> index(dim_);
729  for(int i=0; i<dim_; ++i){
730  index[i] = i;
731  }
732  SortModLeja(theta_,index);
733 
734  //Add roots if neded.
735  ComputeAddedRoots();
736 
737  }
738  }
739 
740  //Function determines whether added roots are needed and adds them if option is turned on.
741  template <class ScalarType, class MV, class OP>
743  {
744  // Store theta (with cols for real and imag parts of Harmonic Ritz Vals)
745  // as one vector of complex numbers to perform arithmetic:
746  std::vector<std::complex<MagnitudeType>> cmplxHRitz (dim_);
747  for(unsigned int i=0; i<cmplxHRitz.size(); ++i){
748  cmplxHRitz[i] = std::complex<MagnitudeType>( theta_(i,0), theta_(i,1) );
749  }
750 
751  // Compute product of factors (pof) to determine added roots:
752  const MagnitudeType one(1.0);
753  std::vector<MagnitudeType> pof (dim_,one);
754  for(int j=0; j<dim_; ++j) {
755  for(int i=0; i<dim_; ++i) {
756  if(i!=j) {
757  pof[j] = std::abs(pof[j]*(one-(cmplxHRitz[j]/cmplxHRitz[i])));
758  }
759  }
760  }
761 
762  // Compute number of extra roots needed:
763  std::vector<int> extra (dim_);
764  int totalExtra = 0;
765  for(int i=0; i<dim_; ++i){
766  extra[i] = ceil((log10(pof[i])-MagnitudeType(4.0))/MagnitudeType(14.0));
767  if(extra[i] > 0){
768  totalExtra += extra[i];
769  }
770  }
771  if (totalExtra){
772  printer_->stream(Warnings) << "Warning: Need to add " << totalExtra << " extra roots." << std::endl;}
773 
774  // If requested to add roots, append them to the theta matrix:
775  if(addRoots_ && totalExtra>0)
776  {
777  theta_.reshape(dim_+totalExtra,2);
778  // Make a matrix copy for perturbed roots:
779  Teuchos::SerialDenseMatrix<OT,MagnitudeType> thetaPert (Teuchos::Copy, theta_, dim_+totalExtra, 2);
780 
781  //Add extra eigenvalues to matrix and perturb for sort:
782  int count = dim_;
783  for(int i=0; i<dim_; ++i){
784  for(int j=0; j< extra[i]; ++j){
785  theta_(count,0) = theta_(i,0);
786  theta_(count,1) = theta_(i,1);
787  thetaPert(count,0) = theta_(i,0)+(j+MCT::one())*MagnitudeType(5e-8);
788  thetaPert(count,1) = theta_(i,1);
789  ++count;
790  }
791  }
792 
793  // Update polynomial degree:
794  dim_ += totalExtra;
795  if (totalExtra){
796  printer_->stream(Warnings) << "New poly degree is: " << dim_ << std::endl;}
797 
798  // Create a new index and sort perturbed roots:
799  std::vector<int> index2(dim_);
800  for(int i=0; i<dim_; ++i){
801  index2[i] = i;
802  }
803  SortModLeja(thetaPert,index2);
804  //Apply sorting to non-perturbed roots:
805  for(int i=0; i<dim_; ++i)
806  {
807  thetaPert(i,0) = theta_(index2[i],0);
808  thetaPert(i,1) = theta_(index2[i],1);
809  }
810  theta_ = thetaPert;
811 
812  }
813  }
814 
815  // Modified Leja sorting function. Takes a serial dense matrix of M harmonic Ritz values and an index
816  // of values from 0 to M. Returns the sorted values and sorted index, similar to Matlab.
817  template <class ScalarType, class MV, class OP>
818  void GmresPolyOp<ScalarType, MV, OP>::SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index) const
819  {
820  //Sort theta values via Modified Leja Ordering:
821 
822  // Set up blank matrices to track sorting:
823  int dimN = index.size();
824  std::vector<int> newIndex(dimN);
828 
829  //Compute all absolute values and find maximum:
830  for(int i = 0; i < dimN; i++){
831  absVal(i) = hypot(thetaN(i,0), thetaN(i,1));
832  }
833  MagnitudeType * maxPointer = std::max_element(absVal.values(), (absVal.values()+dimN));
834  int maxIndex = int (maxPointer- absVal.values());
835 
836  //Put largest abs value first in the list:
837  sorted(0,0) = thetaN(maxIndex,0);
838  sorted(0,1) = thetaN(maxIndex,1);
839  newIndex[0] = index[maxIndex];
840 
841  int j;
842  // If largest value was complex (for real scalar type) put its conjugate in the next slot.
843  if(sorted(0,1)!= SCT::zero() && !SCT::isComplex)
844  {
845  sorted(1,0) = thetaN(maxIndex,0);
846  sorted(1,1) = -thetaN(maxIndex,1);
847  newIndex[1] = index[maxIndex+1];
848  j = 2;
849  }
850  else
851  {
852  j = 1;
853  }
854 
855  //Sort remaining values:
856  MagnitudeType a, b;
857  while( j < dimN )
858  {
859  //For each value, compute (a log of) a product of differences:
860  for(int i = 0; i < dimN; i++)
861  {
862  prod(i) = MCT::one();
863  for(int k = 0; k < j; k++)
864  {
865  a = thetaN(i,0) - sorted(k,0);
866  b = thetaN(i,1) - sorted(k,1);
867  prod(i) = prod(i) + log10(sqrt(a*a + b*b));
868  }
869  }
870 
871  //Value with largest product goes in the next slot:
872  MagnitudeType * maxPointer = std::max_element(prod.values(), (prod.values()+dimN));
873  int maxIndex = int (maxPointer- prod.values());
874  sorted(j,0) = thetaN(maxIndex,0);
875  sorted(j,1) = thetaN(maxIndex,1);
876  newIndex[j] = index[maxIndex];
877 
878  //If it was complex (and scalar type real) put its conjugate in next slot:
879  if(sorted(j,1)!= SCT::zero() && !SCT::isComplex)
880  {
881  j++;
882  sorted(j,0) = thetaN(maxIndex,0);
883  sorted(j,1) = -thetaN(maxIndex,1);
884  newIndex[j] = index[maxIndex+1];
885  }
886  j++;
887  }
888 
889  //Return sorted values and sorted indices:
890  thetaN = sorted;
891  index = newIndex;
892  } //End Modified Leja ordering
893 
894  template <class ScalarType, class MV, class OP>
895  void GmresPolyOp<ScalarType, MV, OP>::ApplyPoly( const MV& x, MV& y ) const
896  {
897  if (dim_) {
898  if (polyType_ == "Arnoldi")
899  ApplyArnoldiPoly(x, y);
900  else if (polyType_ == "Gmres")
901  ApplyGmresPoly(x, y);
902  else if (polyType_ == "Roots")
903  ApplyRootsPoly(x, y);
904  }
905  else {
906  // Just apply the operator in problem_ to x and return y.
907  problem_->applyOp( x, y );
908  }
909  }
910 
911  template <class ScalarType, class MV, class OP>
912  void GmresPolyOp<ScalarType, MV, OP>::ApplyGmresPoly( const MV& x, MV& y ) const
913  {
914  Teuchos::RCP<MV> AX = MVT::CloneCopy(x);
915  Teuchos::RCP<MV> AX2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
916 
917  // Apply left preconditioner.
918  if (!LP_.is_null()) {
919  Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
920  problem_->applyLeftPrec( *AX, *Xtmp ); // Left precondition x into the first vector
921  AX = Xtmp;
922  }
923 
924  {
925 #ifdef BELOS_TEUCHOS_TIME_MONITOR
926  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
927 #endif
928  MVT::MvAddMv(pCoeff_(0,0), *AX, SCT::zero(), y, y); //y= coeff_i(A^ix)
929  }
930  for( int i=1; i < dim_+1; i++)
931  {
932  Teuchos::RCP<MV> X, Y;
933  if ( i%2 )
934  {
935  X = AX;
936  Y = AX2;
937  }
938  else
939  {
940  X = AX2;
941  Y = AX;
942  }
943  problem_->apply(*X, *Y);
944  {
945 #ifdef BELOS_TEUCHOS_TIME_MONITOR
946  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
947 #endif
948  MVT::MvAddMv(pCoeff_(i,0), *Y, SCT::one(), y, y); //y= coeff_i(A^ix) +y
949  }
950  }
951 
952  // Apply right preconditioner.
953  if (!RP_.is_null()) {
954  Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
955  problem_->applyRightPrec( *Ytmp, y );
956  }
957  }
958 
959  template <class ScalarType, class MV, class OP>
960  void GmresPolyOp<ScalarType, MV, OP>::ApplyRootsPoly( const MV& x, MV& y ) const
961  {
962  MVT::MvInit( y, SCT::zero() ); //Zero out y to take the vector with poly applied.
963  Teuchos::RCP<MV> prod = MVT::CloneCopy(x);
964  Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
965  Teuchos::RCP<MV> Xtmp2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
966 
967  // Apply left preconditioner.
968  if (!LP_.is_null()) {
969  problem_->applyLeftPrec( *prod, *Xtmp ); // Left precondition x into the first vector
970  prod = Xtmp;
971  }
972 
973  int i=0;
974  while(i < dim_-1)
975  {
976  if(theta_(i,1)== SCT::zero() || SCT::isComplex) //Real Harmonic Ritz value or complex scalars
977  {
978  {
979 #ifdef BELOS_TEUCHOS_TIME_MONITOR
980  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
981 #endif
982  MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(i,0), *prod, y); //poly = poly + 1/theta_i * prod
983  }
984  problem_->apply(*prod, *Xtmp); // temp = A*prod
985  {
986 #ifdef BELOS_TEUCHOS_TIME_MONITOR
987  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
988 #endif
989  MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/theta_(i,0), *Xtmp, *prod); //prod = prod - 1/theta_i * temp
990  }
991  i++;
992  }
993  else //Current theta is complex and has a conjugate; combine to preserve real arithmetic
994  {
995  MagnitudeType mod = theta_(i,0)*theta_(i,0) + theta_(i,1)*theta_(i,1); //mod = a^2 + b^2
996  problem_->apply(*prod, *Xtmp); // temp = A*prod
997  {
998 #ifdef BELOS_TEUCHOS_TIME_MONITOR
999  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1000 #endif
1001  MVT::MvAddMv(2*theta_(i,0), *prod, -SCT::one(), *Xtmp, *Xtmp); //temp = 2a*prod-temp
1002  MVT::MvAddMv(SCT::one(), y, SCT::one()/mod, *Xtmp, y); //poly = poly + 1/mod*temp
1003  }
1004  if( i < dim_-2 )
1005  {
1006  problem_->apply(*Xtmp, *Xtmp2); // temp2 = A*temp
1007  {
1008 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1009  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1010 #endif
1011  MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/mod, *Xtmp2, *prod); //prod = prod - 1/mod * temp2
1012  }
1013  }
1014  i = i + 2;
1015  }
1016  }
1017  if(theta_(dim_-1,1)== SCT::zero() || SCT::isComplex)
1018  {
1019 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1020  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1021 #endif
1022  MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(dim_-1,0), *prod, y); //poly = poly + 1/theta_i * prod
1023  }
1024 
1025  // Apply right preconditioner.
1026  if (!RP_.is_null()) {
1027  Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
1028  problem_->applyRightPrec( *Ytmp, y );
1029  }
1030  }
1031 
1032  template <class ScalarType, class MV, class OP>
1033  void GmresPolyOp<ScalarType, MV, OP>::ApplyArnoldiPoly( const MV& x, MV& y ) const
1034  {
1035  // Initialize vector storage.
1036  if (V_.is_null()) {
1037  V_ = MVT::Clone( x, dim_ );
1038  if (!LP_.is_null()) {
1039  wL_ = MVT::Clone( y, 1 );
1040  }
1041  if (!RP_.is_null()) {
1042  wR_ = MVT::Clone( y, 1 );
1043  }
1044  }
1045  //
1046  // Apply polynomial to x.
1047  //
1048  int n = MVT::GetNumberVecs( x );
1049  std::vector<int> idxi(1), idxi2, idxj(1);
1050 
1051  // Select vector x[j].
1052  for (int j=0; j<n; ++j) {
1053 
1054  idxi[0] = 0;
1055  idxj[0] = j;
1056  Teuchos::RCP<const MV> x_view = MVT::CloneView( x, idxj );
1057  Teuchos::RCP<MV> y_view = MVT::CloneViewNonConst( y, idxj );
1058  if (!LP_.is_null()) {
1059  Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1060  problem_->applyLeftPrec( *x_view, *v_curr ); // Left precondition x into the first vector of V
1061  } else {
1062  MVT::SetBlock( *x_view, idxi, *V_ ); // Set x as the first vector of V
1063  }
1064 
1065  for (int i=0; i<dim_-1; ++i) {
1066 
1067  // Get views into the current and next vectors
1068  idxi2.resize(i+1);
1069  for (int ii=0; ii<i+1; ++ii) { idxi2[ii] = ii; }
1070  Teuchos::RCP<const MV> v_prev = MVT::CloneView( *V_, idxi2 );
1071  // the tricks below with wR_ and wL_ (potentially set to v_curr and v_next) unfortunately imply that
1072  // v_curr and v_next must be non-const views.
1073  Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1074  idxi[0] = i+1;
1075  Teuchos::RCP<MV> v_next = MVT::CloneViewNonConst( *V_, idxi );
1076 
1077  //---------------------------------------------
1078  // Apply operator to next vector
1079  //---------------------------------------------
1080  // 1) Apply right preconditioner, if we have one.
1081  if (!RP_.is_null()) {
1082  problem_->applyRightPrec( *v_curr, *wR_ );
1083  } else {
1084  wR_ = v_curr;
1085  }
1086  // 2) Check for left preconditioner, if none exists, point at the next vector.
1087  if (LP_.is_null()) {
1088  wL_ = v_next;
1089  }
1090  // 3) Apply operator A.
1091  problem_->applyOp( *wR_, *wL_ );
1092  // 4) Apply left preconditioner, if we have one.
1093  if (!LP_.is_null()) {
1094  problem_->applyLeftPrec( *wL_, *v_next );
1095  }
1096 
1097  // Compute A*v_curr - v_prev*H(1:i,i)
1099  {
1100 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1101  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1102 #endif
1103  MVT::MvTimesMatAddMv( -SCT::one(), *v_prev, h, SCT::one(), *v_next );
1104  }
1105 
1106  // Scale by H(i+1,i)
1107  MVT::MvScale( *v_next, SCT::one()/H_(i+1,i) );
1108  }
1109 
1110  // Compute output y = V*y_./r0_
1111  if (!RP_.is_null()) {
1112  {
1113 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1114  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1115 #endif
1116  MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *wR_ );
1117  }
1118  problem_->applyRightPrec( *wR_, *y_view );
1119  }
1120  else {
1121 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1122  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1123 #endif
1124  MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *y_view );
1125  }
1126  } // (int j=0; j<n; ++j)
1127  } // end Apply()
1128 } // end Belos namespace
1129 
1130 #endif
1131 
1132 // end of file BelosGmresPolyOp.hpp
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:106
ScalarType * values() const
GmresPolyMv(const Teuchos::RCP< const MV > &mv_in)
Collection of types and exceptions used within the Belos solvers.
const GmresPolyMv * CloneView(const std::vector< int > &index) const
Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs ve...
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
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
void MvScale(const std::vector< ScalarType > &alpha)
Scale each element of the i-th vector in *this with alpha[i].
Class which manages the output and verbosity of the Belos solvers.
static void MvDot(const MV &mv, const MV &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
void MvAddMv(const ScalarType alpha, const MultiVec< ScalarType > &A, const ScalarType beta, const MultiVec< ScalarType > &B)
Replace *this with alpha * A + beta * B.
GmresPolyMv * CloneCopy() const
Create a new MultiVec and copy contents of *this into it (deep copy).
void MvScale(const ScalarType alpha)
Scale each element of the vectors in *this with alpha.
ptrdiff_t GetGlobalLength() const
The number of rows in the multivector.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed...
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< MV > getMV()
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static RCP< Time > getNewCounter(const std::string &name)
void ApplyGmresPoly(const MV &x, MV &y) const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
void generateGmresPoly()
This routine takes the matrix, preconditioner, and vectors from the linear problem as well as the par...
Declaration of basic traits for the multivector type.
void POTRS(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
An implementation of StatusTestResNorm using a family of residual norms.
GmresPolyMv * CloneViewNonConst(const std::vector< int > &index)
Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs ve...
GmresPolyMv * CloneCopy(const std::vector< int > &index) const
Creates a new Belos::MultiVec and copies the selected contents of *this into the new multivector (dee...
void MvRandom()
Fill all the vectors in *this with random numbers.
Structure to contain pointers to GmresIteration state variables.
GmresPolyOp(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem_in, const Teuchos::RCP< Teuchos::ParameterList > &params_in)
Basic contstructor.
Belos::StatusTest class for specifying a maximum number of iterations.
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static const double polyTol
Relative residual tolerance for matrix polynomial construction.
Definition: BelosTypes.hpp:304
static void MvTransMv(const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
Class which defines basic traits for the operator type.
int GetNumberVecs() const
The number of vectors (i.e., columns) in the multivector.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
A factory class for generating StatusTestOutput objects.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
ETrans
Whether to apply the (conjugate) transpose of an operator.
Definition: BelosTypes.hpp:81
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
bool isParameter(const std::string &name) const
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index...
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
A Belos::StatusTest class for specifying a maximum number of iterations.
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
void MvTimesMatAddMv(const ScalarType alpha, const MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta)
Update *this with alpha * A * B + beta * (*this).
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
Alternative run-time polymorphic interface for operators.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int curDim
The current dimension of the reduction.
void ApplyPoly(const MV &x, MV &y) const
This routine takes the MV x and applies the polynomial operator phi(OP) to it resulting in the MV y...
void POTRF(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
void factorWithEquilibration(bool flag)
void ApplyArnoldiPoly(const MV &x, MV &y) const
void MvNorm(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm) const
Compute the norm of each vector in *this.
void solveWithTransposeFlag(Teuchos::ETransp trans)
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
GmresPolyOpOrthoFailure(const std::string &what_arg)
void SetBlock(const MultiVec< ScalarType > &A, const std::vector< int > &index)
Copy the vectors in A to a set of vectors in *this.
void ApplyRootsPoly(const MV &x, MV &y) const
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
Teuchos::RCP< const MV > getConstMV() const
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
Belos&#39;s class for applying the GMRES polynomial operator that is used by the hybrid-GMRES linear solv...
Belos concrete class for performing the block GMRES iteration.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
OrdinalType numCols() const
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
NormType
The type of vector norm to compute.
Definition: BelosTypes.hpp:97
void MvInit(const ScalarType alpha)
Replace each element of the vectors in *this with alpha.
void MvDot(const MultiVec< ScalarType > &A, std::vector< ScalarType > &b) const
Compute the dot product of each column of *this with the corresponding column of A.
void MvTransMv(const ScalarType alpha, const MultiVec< ScalarType > &A, Teuchos::SerialDenseMatrix< int, ScalarType > &B) const
Compute a dense matrix B through the matrix-matrix multiply alpha * A^T * (*this).
int setVectors(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B)
Alternative run-time polymorphic interface for operators.
virtual ~GmresPolyOp()
Destructor.
void GEEV(const char &JOBVL, const char &JOBVR, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *WR, MagnitudeType *WI, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
GmresPolyMv * Clone(const int numvecs) const
Create a new MultiVec with numvecs columns.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
bool isType(const std::string &name) const
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
Interface for multivectors used by Belos&#39; linear solvers.
GmresPolyOp(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem_in)
Given no ParameterList, constructor creates no polynomial and only applies the given operator...
static void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of ...
void generateArnoldiPoly()
This routine takes the matrix, preconditioner, and vectors from the linear problem as well as the par...
static void MvScale(MV &mv, const ScalarType alpha)
Scale each element of the vectors in mv with alpha.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
int shape(OrdinalType numRows, OrdinalType numCols)
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...
static void MvPrint(const MV &mv, std::ostream &os)
Print the mv multi-vector to the os output stream.
GmresPolyOpOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorma...
void MvPrint(std::ostream &os) const
Print *this multivector to the os output stream.
OrdinalType stride() const
void Apply(const MultiVec< ScalarType > &x, MultiVec< ScalarType > &y, ETrans=NOTRANS) const
This routine casts the MultiVec to GmresPolyMv to retrieve the MV. Then the above apply method is cal...
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params_in)
Process the passed in parameters.
OrdinalType numRows() const
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
Interface for multivectors used by Belos&#39; linear solvers.
GmresPolyMv(const Teuchos::RCP< MV > &mv_in)
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.

Generated on Fri Aug 14 2020 10:48:33 for Belos by doxygen 1.8.5