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"
61 
65 #include "BelosStatusTestCombo.hpp"
67 
68 #include "BelosOutputManager.hpp"
69 
70 #include "Teuchos_BLAS.hpp"
71 #include "Teuchos_LAPACK.hpp"
72 #include "Teuchos_as.hpp"
73 #include "Teuchos_RCP.hpp"
78 
79 #ifdef BELOS_TEUCHOS_TIME_MONITOR
80  #include "Teuchos_TimeMonitor.hpp"
81 #endif // BELOS_TEUCHOS_TIME_MONITOR
82 
83 namespace Belos {
84 
91  class GmresPolyOpOrthoFailure : public BelosError {public:
92  GmresPolyOpOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
93  {}};
94 
95  // Create a shell class for the MV, inherited off MultiVec<> that will operate with the GmresPolyOp.
96  template <class ScalarType, class MV>
97  class GmresPolyMv : public MultiVec< ScalarType >
98  {
99  public:
100 
101  GmresPolyMv ( const Teuchos::RCP<MV>& mv_in )
102  : mv_(mv_in)
103  {}
105  {
106  mv_ = Teuchos::rcp_const_cast<MV>( mv_in );
107  }
108  Teuchos::RCP<MV> getMV() { return mv_; }
109  Teuchos::RCP<const MV> getConstMV() const { return mv_; }
110 
111  GmresPolyMv * Clone ( const int numvecs ) const
112  {
113  GmresPolyMv * newMV = new GmresPolyMv( MVT::Clone( *mv_, numvecs ) );
114  return newMV;
115  }
116  GmresPolyMv * CloneCopy () const
117  {
118  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneCopy( *mv_ ) );
119  return newMV;
120  }
121  GmresPolyMv * CloneCopy ( const std::vector<int>& index ) const
122  {
123  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneCopy( *mv_, index ) );
124  return newMV;
125  }
126  GmresPolyMv * CloneViewNonConst ( const std::vector<int>& index )
127  {
128  GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneViewNonConst( *mv_, index ) );
129  return newMV;
130  }
131  const GmresPolyMv * CloneView ( const std::vector<int>& index ) const
132  {
133  const GmresPolyMv * newMV = new GmresPolyMv( MVT::CloneView( *mv_, index ) );
134  return newMV;
135  }
136  ptrdiff_t GetGlobalLength () const { return MVT::GetGlobalLength( *mv_ ); }
137  int GetNumberVecs () const { return MVT::GetNumberVecs( *mv_ ); }
138  void MvTimesMatAddMv (const ScalarType alpha,
139  const MultiVec<ScalarType>& A,
140  const Teuchos::SerialDenseMatrix<int,ScalarType>& B, const ScalarType beta)
141  {
142  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
143  MVT::MvTimesMatAddMv( alpha, *(A_in.getConstMV()), B, beta, *mv_ );
144  }
145  void MvAddMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, const ScalarType beta, const MultiVec<ScalarType>& B )
146  {
147  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
148  const GmresPolyMv<ScalarType,MV>& B_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(B);
149  MVT::MvAddMv( alpha, *(A_in.getConstMV()), beta, *(B_in.getConstMV()), *mv_ );
150  }
151  void MvScale ( const ScalarType alpha ) { MVT::MvScale( *mv_, alpha ); }
152  void MvScale ( const std::vector<ScalarType>& alpha ) { MVT::MvScale( *mv_, alpha ); }
153  void MvTransMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, Teuchos::SerialDenseMatrix<int,ScalarType>& B) const
154  {
155  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
156  MVT::MvTransMv( alpha, *(A_in.getConstMV()), *mv_, B );
157  }
158  void MvDot ( const MultiVec<ScalarType>& A, std::vector<ScalarType>& b ) const
159  {
160  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
161  MVT::MvDot( *(A_in.getConstMV()), *mv_, b );
162  }
163  void MvNorm ( std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec, NormType type = TwoNorm ) const
164  {
165  MVT::MvNorm( *mv_, normvec, type );
166  }
167  void SetBlock ( const MultiVec<ScalarType>& A, const std::vector<int>& index )
168  {
169  const GmresPolyMv<ScalarType,MV>& A_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(A);
170  MVT::SetBlock( *(A_in.getConstMV()), index, *mv_ );
171  }
172  void MvRandom () { MVT::MvRandom( *mv_ ); }
173  void MvInit ( const ScalarType alpha ) { MVT::MvInit( *mv_, alpha ); }
174  void MvPrint ( std::ostream& os ) const { MVT::MvPrint( *mv_, os ); }
175 
176  private:
177 
178  typedef MultiVecTraits<ScalarType,MV> MVT;
179 
180  Teuchos::RCP<MV> mv_;
181 
182  };
183 
194  template <class ScalarType, class MV, class OP>
195  class GmresPolyOp : public Operator<ScalarType> {
196  public:
197 
199 
200 
203  const Teuchos::RCP<Teuchos::ParameterList>& params_in
204  )
205  : problem_(problem_in),
206  params_(params_in),
207  LP_(problem_in->getLeftPrec()),
208  RP_(problem_in->getRightPrec())
209  {
210  setParameters( params_ );
211 
212  polyUpdateLabel_ = label_ + ": Hybrid Gmres: Vector Update";
213 #ifdef BELOS_TEUCHOS_TIME_MONITOR
214  timerPolyUpdate_ = Teuchos::TimeMonitor::getNewCounter(polyUpdateLabel_);
215 #endif // BELOS_TEUCHOS_TIME_MONITOR
216 
217  if (polyType_ == "Arnoldi" || polyType_=="Roots")
219  else if (polyType_ == "Gmres")
221  else
222  TEUCHOS_TEST_FOR_EXCEPTION(polyType_!="Arnoldi"&&polyType_!="Gmres"&&polyType_!="Roots",std::invalid_argument,
223  "Belos::GmresPolyOp: \"Polynomial Type\" must be either \"Arnoldi\", \"Gmres\", or \"Roots\".");
224  }
225 
228  : problem_(problem_in)
229  {
230  // If dimension is zero, it will just apply the operator from problem_in in the Apply method.
231  dim_ = 0;
232  }
233 
235  virtual ~GmresPolyOp() {};
237 
239 
240 
242  void setParameters( const Teuchos::RCP<Teuchos::ParameterList>& params_in );
244 
246 
247 
251  void generateArnoldiPoly();
252 
256  void generateGmresPoly();
257 
259 
261 
262 
268  void ApplyPoly ( const MV& x, MV& y ) const;
269  void ApplyArnoldiPoly ( const MV& x, MV& y ) const;
270  void ApplyGmresPoly ( const MV& x, MV& y ) const;
271  void ApplyRootsPoly ( const MV& x, MV& y ) const;
272 
276  void Apply ( const MultiVec<ScalarType>& x, MultiVec<ScalarType>& y, ETrans /* trans */=NOTRANS ) const
277  {
278  const GmresPolyMv<ScalarType,MV>& x_in = dynamic_cast<const GmresPolyMv<ScalarType,MV>&>(x);
279  GmresPolyMv<ScalarType,MV>& y_in = dynamic_cast<GmresPolyMv<ScalarType,MV>&>(y);
280  ApplyPoly( *(x_in.getConstMV()), *(y_in.getMV()) );
281  }
282 
283  int polyDegree() const { return dim_; }
284 
285  private:
286 
287 #ifdef BELOS_TEUCHOS_TIME_MONITOR
288  Teuchos::RCP<Teuchos::Time> timerPolyUpdate_;
289 #endif // BELOS_TEUCHOS_TIME_MONITOR
290  std::string polyUpdateLabel_;
291 
292  typedef int OT; //Ordinal type
293  typedef MultiVecTraits<ScalarType,MV> MVT;
295  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
297 
298  // Default polynomial parameters
299  static constexpr int maxDegree_default_ = 25;
300  static constexpr int verbosity_default_ = Belos::Errors;
301  static constexpr bool randomRHS_default_ = true;
302  static constexpr const char * label_default_ = "Belos";
303  static constexpr const char * polyType_default_ = "Roots";
304  static constexpr const char * orthoType_default_ = "DGKS";
305  static constexpr std::ostream * outputStream_default_ = &std::cout;
306  static constexpr bool damp_default_ = false;
307  static constexpr bool addRoots_default_ = true;
308 
309  // Variables for generating the polynomial
312  Teuchos::RCP<const OP> LP_, RP_;
313 
314  // Output manager.
316  Teuchos::RCP<std::ostream> outputStream_ = Teuchos::rcp(outputStream_default_,false);
317 
318  // Orthogonalization manager.
320 
321  // Current polynomial parameters
322  MagnitudeType polyTol_ = DefaultSolverParameters::polyTol;
323  int maxDegree_ = maxDegree_default_;
324  int verbosity_ = verbosity_default_;
325  bool randomRHS_ = randomRHS_default_;
326  std::string label_ = label_default_;
327  std::string polyType_ = polyType_default_;
328  std::string orthoType_ = orthoType_default_;
329  int dim_ = 0;
330  bool damp_ = damp_default_;
331  bool addRoots_ = addRoots_default_;
332 
333  // Variables for Arnoldi polynomial
334  mutable Teuchos::RCP<MV> V_, wL_, wR_;
337 
338  // Variables for Gmres polynomial;
339  bool autoDeg = false;
341 
342  // Variables for Roots polynomial:
344 
345  // Modified Leja sorting function. Takes a serial dense matrix of M harmonic Ritz values and an index
346  // of values from 0 to M. Returns the sorted values and sorted index, similar to Matlab.
347  void SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index) const ;
348 
349  //Function determines whether added roots are needed and adds them if option is turned on.
350  void ComputeAddedRoots();
351  };
352 
353  template <class ScalarType, class MV, class OP>
355  {
356  // Check which Gmres polynomial to use
357  if (params_in->isParameter("Polynomial Type")) {
358  polyType_ = params_in->get("Polynomial Type", polyType_default_);
359  }
360 
361  // Check for polynomial convergence tolerance
362  if (params_in->isParameter("Polynomial Tolerance")) {
363  if (params_in->isType<MagnitudeType> ("Polynomial Tolerance")) {
364  polyTol_ = params_in->get ("Polynomial Tolerance",
365  static_cast<MagnitudeType> (DefaultSolverParameters::polyTol));
366  }
367  else {
368  polyTol_ = params_in->get ("Polynomial Tolerance", DefaultSolverParameters::polyTol);
369  }
370  }
371 
372  // Check for maximum polynomial degree
373  if (params_in->isParameter("Maximum Degree")) {
374  maxDegree_ = params_in->get("Maximum Degree", maxDegree_default_);
375  }
376 
377  // Check for maximum polynomial degree
378  if (params_in->isParameter("Random RHS")) {
379  randomRHS_ = params_in->get("Random RHS", randomRHS_default_);
380  }
381 
382  // Check for a change in verbosity level
383  if (params_in->isParameter("Verbosity")) {
384  if (Teuchos::isParameterType<int>(*params_in,"Verbosity")) {
385  verbosity_ = params_in->get("Verbosity", verbosity_default_);
386  }
387  else {
388  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params_in,"Verbosity");
389  }
390  }
391 
392  if (params_in->isParameter("Orthogonalization")) {
393  orthoType_ = params_in->get("Orthogonalization",orthoType_default_);
394  }
395 
396  // Check for timer label
397  if (params_in->isParameter("Timer Label")) {
398  label_ = params_in->get("Timer Label", label_default_);
399  }
400 
401  // Output stream
402  if (params_in->isParameter("Output Stream")) {
403  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params_in,"Output Stream");
404  }
405 
406  // Check for damped polynomial
407  if (params_in->isParameter("Damped Poly")) {
408  damp_ = params_in->get("Damped Poly", damp_default_);
409  }
410 
411  // Check for root-adding
412  if (params_in->isParameter("Add Roots")) {
413  addRoots_ = params_in->get("Add Roots", addRoots_default_);
414  }
415  }
416 
417  template <class ScalarType, class MV, class OP>
419  {
420  Teuchos::RCP< MV > V = MVT::Clone( *problem_->getRHS(), maxDegree_+2 );
421 
422  //Make power basis:
423  std::vector<int> index(1,0);
424  Teuchos::RCP< MV > V0 = MVT::CloneViewNonConst(*V, index);
425  if (randomRHS_)
426  MVT::MvRandom( *V0 );
427  else
428  MVT::Assign( *problem_->getRHS(), *V0 );
429 
430  if ( !LP_.is_null() ) {
431  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
432  problem_->applyLeftPrec( *Vtemp, *V0);
433  }
434  if ( damp_ ) {
435  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
436  problem_->apply( *Vtemp, *V0);
437  }
438 
439  for(int i=0; i< maxDegree_+1; i++)
440  {
441  index[0] = i;
442  Teuchos::RCP< const MV > Vi = MVT::CloneView(*V, index);
443  index[0] = i+1;
444  Teuchos::RCP< MV > Vip1 = MVT::CloneViewNonConst(*V, index);
445  problem_->apply( *Vi, *Vip1);
446  }
447 
448  //Consider AV:
449  Teuchos::Range1D range( 1, maxDegree_+1);
450  Teuchos::RCP< const MV > AV = MVT::CloneView( *V, range);
451 
452  //Make lhs (AV)^T(AV)
453  Teuchos::SerialDenseMatrix< OT, ScalarType > AVtransAV( maxDegree_+1, maxDegree_+1);
454  MVT::MvTransMv( SCT::one(), *AV, *AV, AVtransAV);
455  //This process adds pDeg*pDeg + pDeg inner products that aren't in the final count.
456 
458  int infoInt;
459  bool status = true; //Keep adjusting poly deg when true.
460 
461  dim_ = maxDegree_;
463  while( status && dim_ >= 1)
464  {
465  Teuchos::SerialDenseMatrix< OT, ScalarType > lhstemp(Teuchos::Copy, AVtransAV, dim_+1, dim_+1);
466  lapack.POTRF( 'U', dim_+1, lhstemp.values(), lhstemp.stride(), &infoInt);
467 
468  if(autoDeg == false)
469  {
470  status = false;
471  if(infoInt != 0)
472  {
473  std::cout << "BelosGmresPolyOp.hpp: LAPACK POTRF was not successful!!" << std::endl;
474  std::cout << "Error code: " << infoInt << std::endl;
475  }
476  }
477  else
478  {
479  if(infoInt != 0)
480  {//Had bad factor. Reduce poly degree.
481  dim_--;
482  }
483  else
484  {
485  status = false;
486  }
487  }
488  if(status == false)
489  {
490  lhs = lhstemp;
491  }
492  }
493  if(dim_ == 0)
494  {
495  pCoeff_.shape( 1, 1);
496  pCoeff_(0,0) = SCT::one();
497  std::cout << "Poly Degree is zero. No preconditioner created." << std::endl;
498  }
499  else
500  {
501  pCoeff_.shape( dim_+1, 1);
502  //Get correct submatrix of AV:
503  Teuchos::Range1D rangeSub( 1, dim_+1);
504  Teuchos::RCP< const MV > AVsub = MVT::CloneView( *V, rangeSub);
505 
506  //Compute rhs (AV)^T V0
507  MVT::MvTransMv( SCT::one(), *AVsub, *V0, pCoeff_);
508  lapack.POTRS( 'U', dim_+1, 1, lhs.values(), lhs.stride(), pCoeff_.values(), pCoeff_.stride(), &infoInt);
509  if(infoInt != 0)
510  {
511  std::cout << "BelosGmresPolyOp.hpp: LAPACK POTRS was not successful!!" << std::endl;
512  std::cout << "Error code: " << infoInt << std::endl;
513  }
514  }
515  }
516 
517  template <class ScalarType, class MV, class OP>
519  {
520  std::string polyLabel = label_ + ": GmresPolyOp creation";
521 
522  // Create a copy of the linear problem that has a zero initial guess and random RHS.
523  std::vector<int> idx(1,0);
524  Teuchos::RCP<MV> newX = MVT::Clone( *(problem_->getLHS()), 1 );
525  Teuchos::RCP<MV> newB = MVT::Clone( *(problem_->getRHS()), 1 );
526  MVT::MvInit( *newX, SCT::zero() );
527  if (randomRHS_) {
528  MVT::MvRandom( *newB );
529  }
530  else {
531  MVT::Assign( *(MVT::CloneView(*(problem_->getRHS()), idx)), *newB );
532  }
534  Teuchos::rcp( new LinearProblem<ScalarType,MV,OP>( problem_->getOperator(), newX, newB ) );
535  newProblem->setInitResVec( newB );
536  newProblem->setLeftPrec( problem_->getLeftPrec() );
537  newProblem->setRightPrec( problem_->getRightPrec() );
538  newProblem->setLabel(polyLabel);
539  newProblem->setProblem();
540  newProblem->setLSIndex( idx );
541 
542  // Create a parameter list for the GMRES iteration.
543  Teuchos::ParameterList polyList;
544 
545  // Tell the block solver that the block size is one.
546  polyList.set("Num Blocks",maxDegree_);
547  polyList.set("Block Size",1);
548  polyList.set("Keep Hessenberg", true);
549 
550  // Create output manager.
551  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
552 
553  // Create orthogonalization manager if we need to.
554  if (ortho_.is_null()) {
555  params_->set("Orthogonalization", orthoType_);
557  Teuchos::RCP<Teuchos::ParameterList> paramsOrtho; // can be null
558 
559  ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, polyLabel, paramsOrtho);
560  }
561 
562  // Create a simple status test that either reaches the relative residual tolerance or maximum polynomial size.
565 
566  // Implicit residual test, using the native residual to determine if convergence was achieved.
569  convTst->defineScaleForm( convertStringToScaleType("Norm of RHS"), Belos::TwoNorm );
570 
571  // Convergence test that stops the iteration when either are satisfied.
574 
575  // Create Gmres iteration object to perform one cycle of Gmres.
577  gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(newProblem,printer_,polyTest,ortho_,polyList) );
578 
579  // Create the first block in the current Krylov basis (residual).
580  Teuchos::RCP<MV> V_0 = MVT::CloneCopy( *newB );
581  if ( !LP_.is_null() )
582  newProblem->applyLeftPrec( *newB, *V_0 );
583  if ( damp_ )
584  {
585  Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V_0);
586  newProblem->apply( *Vtemp, *V_0 );
587  }
588 
589  // Get a matrix to hold the orthonormalization coefficients.
590  r0_.resize(1);
591 
592  // Orthonormalize the new V_0
593  int rank = ortho_->normalize( *V_0, Teuchos::rcpFromRef(r0_) );
595  "Belos::GmresPolyOp::generateArnoldiPoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
596 
597  // Set the new state and initialize the solver.
599  newstate.V = V_0;
600  newstate.z = Teuchos::rcpFromRef( r0_);
601  newstate.curDim = 0;
602  gmres_iter->initializeGmres(newstate);
603 
604  // Perform Gmres iteration
605  try {
606  gmres_iter->iterate();
607  }
608  catch (GmresIterationOrthoFailure& e) {
609  // Try to recover the most recent least-squares solution
610  gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
611  }
612  catch (std::exception& e) {
613  using std::endl;
614  printer_->stream(Errors) << "Error! Caught exception in BlockGmresIter::iterate() at iteration "
615  << gmres_iter->getNumIters() << endl << e.what () << endl;
616  throw;
617  }
618 
619  // Get the solution for this polynomial, use in comparison below
620  Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
621 
622  // Record polynomial info, get current GMRES state
623  GmresIterationState<ScalarType,MV> gmresState = gmres_iter->getState();
624 
625  // If the polynomial has no dimension, the tolerance is too low, return false
626  dim_ = gmresState.curDim;
627  if (dim_ == 0) {
628  return;
629  }
630  if(polyType_ == "Arnoldi"){
631  // Make a view and then copy the RHS of the least squares problem.
632  //
633  y_ = Teuchos::SerialDenseMatrix<OT,ScalarType>( Teuchos::Copy, *gmresState.z, dim_, 1 );
634  H_ = *gmresState.H;
635 
636  //
637  // Solve the least squares problem.
638  //
641  Teuchos::NON_UNIT_DIAG, dim_, 1, SCT::one(),
642  gmresState.R->values(), gmresState.R->stride(),
643  y_.values(), y_.stride() );
644  }
645  else{ //Generate Roots Poly
646  //Find Harmonic Ritz Values to use as polynomial roots:
647 
648  //Copy of square H used to find poly roots:
649  H_ = Teuchos::SerialDenseMatrix<OT,ScalarType>(Teuchos::Copy, *gmresState.H, dim_, dim_);
650  //Zero out below subdiagonal of H:
651  for(int i=0; i <= dim_-3; i++) {
652  for(int k=i+2; k <= dim_-1; k++) {
653  H_(k,i) = SCT::zero();
654  }
655  }
656  //Extra copy of H because equilibrate changes the matrix:
657  Teuchos::SerialDenseMatrix<OT,ScalarType> Htemp (Teuchos::Copy, H_, dim_, dim_);
658 
659  //View the m+1,m element and last col of H:
660  ScalarType Hlast = (*gmresState.H)(dim_,dim_-1);
661  Teuchos::SerialDenseMatrix<OT,ScalarType> HlastCol (Teuchos::View, H_, dim_, 1, 0, dim_-1);
662 
663  //Set up linear system for H^{-*}e_m:
665  E.putScalar(SCT::zero());
666  E(dim_-1,0) = SCT::one();
667 
669  HSolver.setMatrix( Teuchos::rcpFromRef(Htemp));
671  HSolver.setVectors( Teuchos::rcpFromRef(F), Teuchos::rcpFromRef(E));
672  HSolver.factorWithEquilibration( true );
673 
674  //Factor matrix and solve for F = H^{-*}e_m:
675  int info = 0;
676  info = HSolver.factor();
677  if(info != 0){
678  std::cout << "Hsolver factor: info = " << info << std::endl;
679  }
680  info = HSolver.solve();
681  if(info != 0){
682  std::cout << "Hsolver solve : info = " << info << std::endl;
683  }
684 
685  //Scale F and adjust H for Harmonic Ritz value eigenproblem:
686  F.scale(Hlast*Hlast);
687  HlastCol += F;
688 
689  //Set up for eigenvalue problem to get Harmonic Ritz Values:
691  theta_.shape(dim_,2);//1st col for real part, 2nd col for imaginary
692 
693  const int ldv = 1;
694  ScalarType* vlr = 0;
695 
696  // Size of workspace and workspace for DGEEV
697  int lwork = -1;
698  std::vector<ScalarType> work(1);
699  std::vector<MagnitudeType> rwork(2*dim_);
700 
701  //Find workspace size for DGEEV:
702  lapack.GEEV('N','N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
703  lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
704  work.resize( lwork );
705  // Solve for Harmonic Ritz Values:
706  lapack.GEEV('N','N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
707 
708  if(info != 0){
709  std::cout << "GEEV solve : info = " << info << std::endl;
710  }
711 
712  // Set index for sort function and sort Harmonic Ritz Values:
713  std::vector<int> index(dim_);
714  for(int i=0; i<dim_; ++i){
715  index[i] = i;
716  }
717  SortModLeja(theta_,index);
718 
719  //Add roots if neded.
720  ComputeAddedRoots();
721 
722  }
723  }
724 
725  //Function determines whether added roots are needed and adds them if option is turned on.
726  template <class ScalarType, class MV, class OP>
728  {
729  // Store theta (with cols for real and imag parts of Harmonic Ritz Vals)
730  // as one vector of complex numbers to perform arithmetic:
731  std::vector<std::complex<MagnitudeType>> cmplxHRitz (dim_);
732  for(unsigned int i=0; i<cmplxHRitz.size(); ++i){
733  cmplxHRitz[i] = std::complex<MagnitudeType>( theta_(i,0), theta_(i,1) );
734  }
735 
736  // Compute product of factors (pof) to determine added roots:
737  const MagnitudeType one(1.0);
738  std::vector<MagnitudeType> pof (dim_,one);
739  for(int j=0; j<dim_; ++j) {
740  for(int i=0; i<dim_; ++i) {
741  if(i!=j) {
742  pof[j] = std::abs(pof[j]*(one-(cmplxHRitz[j]/cmplxHRitz[i])));
743  }
744  }
745  }
746 
747  // Compute number of extra roots needed:
748  std::vector<int> extra (dim_);
749  int totalExtra = 0;
750  for(int i=0; i<dim_; ++i){
751  if (pof[i] > MCT::zero())
752  extra[i] = ceil((log10(pof[i])-MagnitudeType(4.0))/MagnitudeType(14.0));
753  else
754  extra[i] = 0;
755  if(extra[i] > 0){
756  totalExtra += extra[i];
757  }
758  }
759  if (totalExtra){
760  printer_->stream(Warnings) << "Warning: Need to add " << totalExtra << " extra roots." << std::endl;}
761 
762  // If requested to add roots, append them to the theta matrix:
763  if(addRoots_ && totalExtra>0)
764  {
765  theta_.reshape(dim_+totalExtra,2);
766  // Make a matrix copy for perturbed roots:
767  Teuchos::SerialDenseMatrix<OT,MagnitudeType> thetaPert (Teuchos::Copy, theta_, dim_+totalExtra, 2);
768 
769  //Add extra eigenvalues to matrix and perturb for sort:
770  int count = dim_;
771  for(int i=0; i<dim_; ++i){
772  for(int j=0; j< extra[i]; ++j){
773  theta_(count,0) = theta_(i,0);
774  theta_(count,1) = theta_(i,1);
775  thetaPert(count,0) = theta_(i,0)+(j+MCT::one())*MagnitudeType(5e-8);
776  thetaPert(count,1) = theta_(i,1);
777  ++count;
778  }
779  }
780 
781  // Update polynomial degree:
782  dim_ += totalExtra;
783  if (totalExtra){
784  printer_->stream(Warnings) << "New poly degree is: " << dim_ << std::endl;}
785 
786  // Create a new index and sort perturbed roots:
787  std::vector<int> index2(dim_);
788  for(int i=0; i<dim_; ++i){
789  index2[i] = i;
790  }
791  SortModLeja(thetaPert,index2);
792  //Apply sorting to non-perturbed roots:
793  for(int i=0; i<dim_; ++i)
794  {
795  thetaPert(i,0) = theta_(index2[i],0);
796  thetaPert(i,1) = theta_(index2[i],1);
797  }
798  theta_ = thetaPert;
799 
800  }
801  }
802 
803  // Modified Leja sorting function. Takes a serial dense matrix of M harmonic Ritz values and an index
804  // of values from 0 to M. Returns the sorted values and sorted index, similar to Matlab.
805  template <class ScalarType, class MV, class OP>
806  void GmresPolyOp<ScalarType, MV, OP>::SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index) const
807  {
808  //Sort theta values via Modified Leja Ordering:
809 
810  // Set up blank matrices to track sorting:
811  int dimN = index.size();
812  std::vector<int> newIndex(dimN);
816 
817  //Compute all absolute values and find maximum:
818  for(int i = 0; i < dimN; i++){
819  absVal(i) = hypot(thetaN(i,0), thetaN(i,1));
820  }
821  MagnitudeType * maxPointer = std::max_element(absVal.values(), (absVal.values()+dimN));
822  int maxIndex = int (maxPointer- absVal.values());
823 
824  //Put largest abs value first in the list:
825  sorted(0,0) = thetaN(maxIndex,0);
826  sorted(0,1) = thetaN(maxIndex,1);
827  newIndex[0] = index[maxIndex];
828 
829  int j;
830  // If largest value was complex (for real scalar type) put its conjugate in the next slot.
831  if(sorted(0,1)!= SCT::zero() && !SCT::isComplex)
832  {
833  sorted(1,0) = thetaN(maxIndex,0);
834  sorted(1,1) = -thetaN(maxIndex,1);
835  newIndex[1] = index[maxIndex+1];
836  j = 2;
837  }
838  else
839  {
840  j = 1;
841  }
842 
843  //Sort remaining values:
844  MagnitudeType a, b;
845  while( j < dimN )
846  {
847  //For each value, compute (a log of) a product of differences:
848  for(int i = 0; i < dimN; i++)
849  {
850  prod(i) = MCT::one();
851  for(int k = 0; k < j; k++)
852  {
853  a = thetaN(i,0) - sorted(k,0);
854  b = thetaN(i,1) - sorted(k,1);
855  if (a*a + b*b > MCT::zero())
856  prod(i) = prod(i) + log10(hypot(a,b));
857  else {
858  prod(i) = -std::numeric_limits<MagnitudeType>::infinity();
859  break;
860  }
861  }
862  }
863 
864  //Value with largest product goes in the next slot:
865  MagnitudeType * maxPointer = std::max_element(prod.values(), (prod.values()+dimN));
866  int maxIndex = int (maxPointer- prod.values());
867  sorted(j,0) = thetaN(maxIndex,0);
868  sorted(j,1) = thetaN(maxIndex,1);
869  newIndex[j] = index[maxIndex];
870 
871  //If it was complex (and scalar type real) put its conjugate in next slot:
872  if(sorted(j,1)!= SCT::zero() && !SCT::isComplex)
873  {
874  j++;
875  sorted(j,0) = thetaN(maxIndex,0);
876  sorted(j,1) = -thetaN(maxIndex,1);
877  newIndex[j] = index[maxIndex+1];
878  }
879  j++;
880  }
881 
882  //Return sorted values and sorted indices:
883  thetaN = sorted;
884  index = newIndex;
885  } //End Modified Leja ordering
886 
887  template <class ScalarType, class MV, class OP>
888  void GmresPolyOp<ScalarType, MV, OP>::ApplyPoly( const MV& x, MV& y ) const
889  {
890  if (dim_) {
891  if (polyType_ == "Arnoldi")
892  ApplyArnoldiPoly(x, y);
893  else if (polyType_ == "Gmres")
894  ApplyGmresPoly(x, y);
895  else if (polyType_ == "Roots")
896  ApplyRootsPoly(x, y);
897  }
898  else {
899  // Just apply the operator in problem_ to x and return y.
900  problem_->applyOp( x, y );
901  }
902  }
903 
904  template <class ScalarType, class MV, class OP>
905  void GmresPolyOp<ScalarType, MV, OP>::ApplyGmresPoly( const MV& x, MV& y ) const
906  {
907  Teuchos::RCP<MV> AX = MVT::CloneCopy(x);
908  Teuchos::RCP<MV> AX2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
909 
910  // Apply left preconditioner.
911  if (!LP_.is_null()) {
912  Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
913  problem_->applyLeftPrec( *AX, *Xtmp ); // Left precondition x into the first vector
914  AX = Xtmp;
915  }
916 
917  {
918 #ifdef BELOS_TEUCHOS_TIME_MONITOR
919  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
920 #endif
921  MVT::MvAddMv(pCoeff_(0,0), *AX, SCT::zero(), y, y); //y= coeff_i(A^ix)
922  }
923  for( int i=1; i < dim_+1; i++)
924  {
925  Teuchos::RCP<MV> X, Y;
926  if ( i%2 )
927  {
928  X = AX;
929  Y = AX2;
930  }
931  else
932  {
933  X = AX2;
934  Y = AX;
935  }
936  problem_->apply(*X, *Y);
937  {
938 #ifdef BELOS_TEUCHOS_TIME_MONITOR
939  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
940 #endif
941  MVT::MvAddMv(pCoeff_(i,0), *Y, SCT::one(), y, y); //y= coeff_i(A^ix) +y
942  }
943  }
944 
945  // Apply right preconditioner.
946  if (!RP_.is_null()) {
947  Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
948  problem_->applyRightPrec( *Ytmp, y );
949  }
950  }
951 
952  template <class ScalarType, class MV, class OP>
953  void GmresPolyOp<ScalarType, MV, OP>::ApplyRootsPoly( const MV& x, MV& y ) const
954  {
955  MVT::MvInit( y, SCT::zero() ); //Zero out y to take the vector with poly applied.
956  Teuchos::RCP<MV> prod = MVT::CloneCopy(x);
957  Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
958  Teuchos::RCP<MV> Xtmp2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
959 
960  // Apply left preconditioner.
961  if (!LP_.is_null()) {
962  problem_->applyLeftPrec( *prod, *Xtmp ); // Left precondition x into the first vector
963  prod = Xtmp;
964  }
965 
966  int i=0;
967  while(i < dim_-1)
968  {
969  if(theta_(i,1)== SCT::zero() || SCT::isComplex) //Real Harmonic Ritz value or complex scalars
970  {
971  {
972 #ifdef BELOS_TEUCHOS_TIME_MONITOR
973  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
974 #endif
975  MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(i,0), *prod, y); //poly = poly + 1/theta_i * prod
976  }
977  problem_->apply(*prod, *Xtmp); // temp = A*prod
978  {
979 #ifdef BELOS_TEUCHOS_TIME_MONITOR
980  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
981 #endif
982  MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/theta_(i,0), *Xtmp, *prod); //prod = prod - 1/theta_i * temp
983  }
984  i++;
985  }
986  else //Current theta is complex and has a conjugate; combine to preserve real arithmetic
987  {
988  MagnitudeType mod = theta_(i,0)*theta_(i,0) + theta_(i,1)*theta_(i,1); //mod = a^2 + b^2
989  problem_->apply(*prod, *Xtmp); // temp = A*prod
990  {
991 #ifdef BELOS_TEUCHOS_TIME_MONITOR
992  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
993 #endif
994  MVT::MvAddMv(2*theta_(i,0), *prod, -SCT::one(), *Xtmp, *Xtmp); //temp = 2a*prod-temp
995  MVT::MvAddMv(SCT::one(), y, SCT::one()/mod, *Xtmp, y); //poly = poly + 1/mod*temp
996  }
997  if( i < dim_-2 )
998  {
999  problem_->apply(*Xtmp, *Xtmp2); // temp2 = A*temp
1000  {
1001 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1002  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1003 #endif
1004  MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/mod, *Xtmp2, *prod); //prod = prod - 1/mod * temp2
1005  }
1006  }
1007  i = i + 2;
1008  }
1009  }
1010  if(theta_(dim_-1,1)== SCT::zero() || SCT::isComplex)
1011  {
1012 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1013  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1014 #endif
1015  MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(dim_-1,0), *prod, y); //poly = poly + 1/theta_i * prod
1016  }
1017 
1018  // Apply right preconditioner.
1019  if (!RP_.is_null()) {
1020  Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
1021  problem_->applyRightPrec( *Ytmp, y );
1022  }
1023  }
1024 
1025  template <class ScalarType, class MV, class OP>
1026  void GmresPolyOp<ScalarType, MV, OP>::ApplyArnoldiPoly( const MV& x, MV& y ) const
1027  {
1028  // Initialize vector storage.
1029  if (V_.is_null()) {
1030  V_ = MVT::Clone( x, dim_ );
1031  if (!LP_.is_null()) {
1032  wL_ = MVT::Clone( y, 1 );
1033  }
1034  if (!RP_.is_null()) {
1035  wR_ = MVT::Clone( y, 1 );
1036  }
1037  }
1038  //
1039  // Apply polynomial to x.
1040  //
1041  int n = MVT::GetNumberVecs( x );
1042  std::vector<int> idxi(1), idxi2, idxj(1);
1043 
1044  // Select vector x[j].
1045  for (int j=0; j<n; ++j) {
1046 
1047  idxi[0] = 0;
1048  idxj[0] = j;
1049  Teuchos::RCP<const MV> x_view = MVT::CloneView( x, idxj );
1050  Teuchos::RCP<MV> y_view = MVT::CloneViewNonConst( y, idxj );
1051  if (!LP_.is_null()) {
1052  Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1053  problem_->applyLeftPrec( *x_view, *v_curr ); // Left precondition x into the first vector of V
1054  } else {
1055  MVT::SetBlock( *x_view, idxi, *V_ ); // Set x as the first vector of V
1056  }
1057 
1058  for (int i=0; i<dim_-1; ++i) {
1059 
1060  // Get views into the current and next vectors
1061  idxi2.resize(i+1);
1062  for (int ii=0; ii<i+1; ++ii) { idxi2[ii] = ii; }
1063  Teuchos::RCP<const MV> v_prev = MVT::CloneView( *V_, idxi2 );
1064  // the tricks below with wR_ and wL_ (potentially set to v_curr and v_next) unfortunately imply that
1065  // v_curr and v_next must be non-const views.
1066  Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1067  idxi[0] = i+1;
1068  Teuchos::RCP<MV> v_next = MVT::CloneViewNonConst( *V_, idxi );
1069 
1070  //---------------------------------------------
1071  // Apply operator to next vector
1072  //---------------------------------------------
1073  // 1) Apply right preconditioner, if we have one.
1074  if (!RP_.is_null()) {
1075  problem_->applyRightPrec( *v_curr, *wR_ );
1076  } else {
1077  wR_ = v_curr;
1078  }
1079  // 2) Check for left preconditioner, if none exists, point at the next vector.
1080  if (LP_.is_null()) {
1081  wL_ = v_next;
1082  }
1083  // 3) Apply operator A.
1084  problem_->applyOp( *wR_, *wL_ );
1085  // 4) Apply left preconditioner, if we have one.
1086  if (!LP_.is_null()) {
1087  problem_->applyLeftPrec( *wL_, *v_next );
1088  }
1089 
1090  // Compute A*v_curr - v_prev*H(1:i,i)
1092  {
1093 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1094  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1095 #endif
1096  MVT::MvTimesMatAddMv( -SCT::one(), *v_prev, h, SCT::one(), *v_next );
1097  }
1098 
1099  // Scale by H(i+1,i)
1100  MVT::MvScale( *v_next, SCT::one()/H_(i+1,i) );
1101  }
1102 
1103  // Compute output y = V*y_./r0_
1104  if (!RP_.is_null()) {
1105  {
1106 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1107  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1108 #endif
1109  MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *wR_ );
1110  }
1111  problem_->applyRightPrec( *wR_, *y_view );
1112  }
1113  else {
1114 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1115  Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1116 #endif
1117  MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *y_view );
1118  }
1119  } // (int j=0; j<n; ++j)
1120  } // end Apply()
1121 } // end Belos namespace
1122 
1123 #endif
1124 
1125 // 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:296
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.
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...
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.
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
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.
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
Return an instance of the specified MatOrthoManager subclass.

Generated on Thu Nov 21 2024 09:28:19 for Belos by doxygen 1.8.5