Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBlockGmresIter.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_BLOCK_GMRES_ITER_HPP
11 #define BELOS_BLOCK_GMRES_ITER_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosGmresIteration.hpp"
20 
21 #include "BelosLinearProblem.hpp"
22 #include "BelosMatOrthoManager.hpp"
23 #include "BelosOutputManager.hpp"
24 #include "BelosStatusTest.hpp"
25 #include "BelosOperatorTraits.hpp"
26 #include "BelosMultiVecTraits.hpp"
27 
28 #include "Teuchos_BLAS.hpp"
29 #include "Teuchos_LAPACK.hpp"
32 #include "Teuchos_ScalarTraits.hpp"
34 #include "Teuchos_TimeMonitor.hpp"
35 
49 namespace Belos {
50 
51 template<class ScalarType, class MV, class OP>
52 class BlockGmresIter : virtual public GmresIteration<ScalarType,MV,OP> {
53 
54  public:
55 
56  //
57  // Convenience typedefs
58  //
63 
65 
66 
77  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
80  Teuchos::ParameterList &params );
81 
83  virtual ~BlockGmresIter() {};
85 
86 
88 
89 
111  void iterate();
112 
135 
139  void initialize()
140  {
142  initializeGmres(empty);
143  }
144 
154  state.curDim = curDim_;
155  state.V = V_;
156  state.H = H_;
157  state.R = R_;
158  state.z = z_;
159  return state;
160  }
161 
163 
164 
166 
167 
169  int getNumIters() const { return iter_; }
170 
172  void resetNumIters( int iter = 0 ) { iter_ = iter; }
173 
176  Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
177 
179 
185 
187 
190  void updateLSQR( int dim = -1 );
191 
193  int getCurSubspaceDim() const {
194  if (!initialized_) return 0;
195  return curDim_;
196  };
197 
199  int getMaxSubspaceDim() const { return blockSize_*numBlocks_; }
200 
202 
203 
205 
206 
208  const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
209 
211  int getBlockSize() const { return blockSize_; }
212 
214  void setBlockSize(int blockSize) { setSize( blockSize, numBlocks_ ); }
215 
217  int getNumBlocks() const { return numBlocks_; }
218 
220  void setNumBlocks(int numBlocks) { setSize( blockSize_, numBlocks ); }
221 
228  void setSize(int blockSize, int numBlocks);
229 
231  bool isInitialized() { return initialized_; }
232 
234 
235  private:
236 
237  //
238  // Internal structs
239  //
240  struct CheckList {
241  bool checkV;
242  bool checkArn;
243  CheckList() : checkV(false), checkArn(false) {};
244  };
245  //
246  // Internal methods
247  //
249  std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
250 
252  void setStateSize();
253 
254  //
255  // Classes inputed through constructor that define the linear problem to be solved.
256  //
261 
262  //
263  // Algorithmic parameters
264  //
265  // blockSize_ is the solver block size.
266  // It controls the number of vectors added to the basis on each iteration.
267  int blockSize_;
268  // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
269  int numBlocks_;
270 
271  // Storage for QR factorization of the least squares system.
274 
275  //
276  // Current solver state
277  //
278  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
279  // is capable of running; _initialize is controlled by the initialize() member method
280  // For the implications of the state of initialized_, please see documentation for initialize()
281  bool initialized_;
282 
283  // stateStorageInitialized_ specifies that the state storage has be initialized to the current
284  // blockSize_ and numBlocks_. This initialization may be postponed if the linear problem was
285  // generated without the right-hand side or solution vectors.
286  bool stateStorageInitialized_;
287 
288  // keepHessenberg_ specifies that the iteration must keep the Hessenberg matrix formed via the
289  // Arnoldi factorization and the upper triangular matrix that is the Hessenberg matrix reduced via
290  // QR factorization separate.
291  bool keepHessenberg_;
292 
293  // initHessenberg_ specifies that the iteration should reinitialize the Hessenberg matrix by zeroing
294  // out all entries before an iteration is started.
295  bool initHessenberg_;
296 
297  // Current subspace dimension, and number of iterations performed.
298  int curDim_, iter_;
299 
300  //
301  // State Storage
302  //
303  Teuchos::RCP<MV> V_;
304  //
305  // Projected matrices
306  // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
307  //
309  //
310  // QR decomposition of Projected matrices for solving the least squares system HY = B.
311  // R_: Upper triangular reduction of H
312  // z_: Q applied to right-hand side of the least squares system
315 };
316 
318  // Constructor.
319  template<class ScalarType, class MV, class OP>
321  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
324  Teuchos::ParameterList &params ):
325  lp_(problem),
326  om_(printer),
327  stest_(tester),
328  ortho_(ortho),
329  blockSize_(0),
330  numBlocks_(0),
331  initialized_(false),
332  stateStorageInitialized_(false),
333  keepHessenberg_(false),
334  initHessenberg_(false),
335  curDim_(0),
336  iter_(0)
337  {
338  // Find out whether we are saving the Hessenberg matrix.
339  if ( om_->isVerbosity( Debug ) )
340  keepHessenberg_ = true;
341  else
342  keepHessenberg_ = params.get("Keep Hessenberg", false);
343 
344  // Find out whether we are initializing the Hessenberg matrix.
345  initHessenberg_ = params.get("Initialize Hessenberg", false);
346 
347  // Get the maximum number of blocks allowed for this Krylov subspace
348  TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter("Num Blocks"), std::invalid_argument,
349  "Belos::BlockGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
350  int nb = Teuchos::getParameter<int>(params, "Num Blocks");
351 
352  // Set the block size and allocate data
353  int bs = params.get("Block Size", 1);
354  setSize( bs, nb );
355  }
356 
358  // Set the block size and make necessary adjustments.
359  template <class ScalarType, class MV, class OP>
360  void BlockGmresIter<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
361  {
362  // This routine only allocates space; it doesn't not perform any computation
363  // any change in size will invalidate the state of the solver.
364 
365  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Belos::BlockGmresIter::setSize was passed a non-positive argument.");
366  if (blockSize == blockSize_ && numBlocks == numBlocks_) {
367  // do nothing
368  return;
369  }
370 
371  if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
372  stateStorageInitialized_ = false;
373 
374  blockSize_ = blockSize;
375  numBlocks_ = numBlocks;
376 
377  initialized_ = false;
378  curDim_ = 0;
379 
380  // Use the current blockSize_ and numBlocks_ to initialize the state storage.
381  setStateSize();
382 
383  }
384 
386  // Setup the state storage.
387  template <class ScalarType, class MV, class OP>
389  {
390  if (!stateStorageInitialized_) {
391 
392  // Check if there is any multivector to clone from.
393  Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
394  Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
395  if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
396  stateStorageInitialized_ = false;
397  return;
398  }
399  else {
400 
402  // blockSize*numBlocks dependent
403  //
404  int newsd = blockSize_*(numBlocks_+1);
405 
406  if (blockSize_==1) {
407  cs.resize( newsd );
408  sn.resize( newsd );
409  }
410  else {
411  beta.resize( newsd );
412  }
413 
414  // Initialize the state storage
415  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_*static_cast<ptrdiff_t>(numBlocks_) > MVT::GetGlobalLength(*rhsMV),std::invalid_argument,
416  "Belos::BlockGmresIter::setStateSize(): Cannot generate a Krylov basis with dimension larger the operator!");
417 
418  // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
419  if (V_ == Teuchos::null) {
420  // Get the multivector that is not null.
421  Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
422  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
423  "Belos::BlockGmresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
424  V_ = MVT::Clone( *tmp, newsd );
425  }
426  else {
427  // Generate V_ by cloning itself ONLY if more space is needed.
428  if (MVT::GetNumberVecs(*V_) < newsd) {
429  Teuchos::RCP<const MV> tmp = V_;
430  V_ = MVT::Clone( *tmp, newsd );
431  }
432  }
433 
434  // Generate R_ only if it doesn't exist, otherwise resize it.
435  if (R_ == Teuchos::null) {
437  }
438  if (initHessenberg_) {
439  R_->shape( newsd, newsd-blockSize_ );
440  }
441  else {
442  if (R_->numRows() < newsd || R_->numCols() < newsd-blockSize_) {
443  R_->shapeUninitialized( newsd, newsd-blockSize_ );
444  }
445  }
446 
447  // Generate H_ only if it doesn't exist, and we are keeping the upper Hessenberg matrix.
448  if (keepHessenberg_) {
449  if (H_ == Teuchos::null) {
451  }
452  if (initHessenberg_) {
453  H_->shape( newsd, newsd-blockSize_ );
454  }
455  else {
456  if (H_->numRows() < newsd || H_->numCols() < newsd-blockSize_) {
457  H_->shapeUninitialized( newsd, newsd-blockSize_ );
458  }
459  }
460  }
461  else {
462  // Point H_ and R_ at the same object.
463  H_ = R_;
464  }
465 
466  // Generate z_ only if it doesn't exist, otherwise resize it.
467  if (z_ == Teuchos::null) {
469  }
470  if (z_-> numRows() < newsd || z_->numCols() < blockSize_) {
471  z_->shapeUninitialized( newsd, blockSize_ );
472  }
473 
474  // State storage has now been initialized.
475  stateStorageInitialized_ = true;
476  }
477  }
478  }
479 
481  // Get the current update from this subspace.
482  template <class ScalarType, class MV, class OP>
484  {
485  //
486  // If this is the first iteration of the Arnoldi factorization,
487  // there is no update, so return Teuchos::null.
488  //
489  Teuchos::RCP<MV> currentUpdate = Teuchos::null;
490  if (curDim_==0) {
491  return currentUpdate;
492  } else {
493  const ScalarType one = SCT::one();
494  const ScalarType zero = SCT::zero();
496  currentUpdate = MVT::Clone( *V_, blockSize_ );
497  //
498  // Make a view and then copy the RHS of the least squares problem. DON'T OVERWRITE IT!
499  //
500  Teuchos::SerialDenseMatrix<int,ScalarType> y( Teuchos::Copy, *z_, curDim_, blockSize_ );
501  //
502  // Solve the least squares problem.
503  //
505  Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,
506  R_->values(), R_->stride(), y.values(), y.stride() );
507  //
508  // Compute the current update.
509  //
510  std::vector<int> index(curDim_);
511  for ( int i=0; i<curDim_; i++ ) {
512  index[i] = i;
513  }
514  Teuchos::RCP<const MV> Vjp1 = MVT::CloneView( *V_, index );
515  MVT::MvTimesMatAddMv( one, *Vjp1, y, zero, *currentUpdate );
516  }
517  return currentUpdate;
518  }
519 
520 
522  // Get the native residuals stored in this iteration.
523  // Note: No residual std::vector will be returned by Gmres.
524  template <class ScalarType, class MV, class OP>
526  {
527  //
528  // NOTE: Make sure the incoming std::vector is the correct size!
529  //
530  if ( norms && (int)norms->size() < blockSize_ )
531  norms->resize( blockSize_ );
532 
533  if (norms) {
535  for (int j=0; j<blockSize_; j++) {
536  (*norms)[j] = blas.NRM2( blockSize_, &(*z_)(curDim_, j), 1);
537  }
538  }
539  return Teuchos::null;
540  }
541 
542 
543 
545  // Initialize this iteration object
546  template <class ScalarType, class MV, class OP>
548  {
549  // Initialize the state storage if it isn't already.
550  if (!stateStorageInitialized_)
551  setStateSize();
552 
553  TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
554  "Belos::BlockGmresIter::initialize(): Cannot initialize state storage!");
555 
556  const ScalarType zero = SCT::zero(), one = SCT::one();
557 
558  // NOTE: In BlockGmresIter, V and Z are required!!!
559  // inconsitent multivectors widths and lengths will not be tolerated, and
560  // will be treated with exceptions.
561  //
562  std::string errstr("Belos::BlockGmresIter::initialize(): Specified multivectors must have a consistent length and width.");
563 
564  if (newstate.V != Teuchos::null && newstate.z != Teuchos::null) {
565 
566  // initialize V_,z_, and curDim_
567 
568  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_),
569  std::invalid_argument, errstr );
570  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < blockSize_,
571  std::invalid_argument, errstr );
572  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*(numBlocks_+1),
573  std::invalid_argument, errstr );
574 
575  curDim_ = newstate.curDim;
576  int lclDim = MVT::GetNumberVecs(*newstate.V);
577 
578  // check size of Z
579  TEUCHOS_TEST_FOR_EXCEPTION(newstate.z->numRows() < curDim_ || newstate.z->numCols() < blockSize_, std::invalid_argument, errstr);
580 
581 
582  // copy basis vectors from newstate into V
583  if (newstate.V != V_) {
584  // only copy over the first block and print a warning.
585  if (curDim_ == 0 && lclDim > blockSize_) {
586  om_->stream(Warnings) << "Belos::BlockGmresIter::initialize(): the solver was initialized with a kernel of " << lclDim << std::endl
587  << "The block size however is only " << blockSize_ << std::endl
588  << "The last " << lclDim - blockSize_ << " vectors will be discarded." << std::endl;
589  }
590  std::vector<int> nevind(curDim_+blockSize_);
591  for (int i=0; i<curDim_+blockSize_; i++) nevind[i] = i;
592  Teuchos::RCP<const MV> newV = MVT::CloneView( *newstate.V, nevind );
593  Teuchos::RCP<MV> lclV = MVT::CloneViewNonConst( *V_, nevind );
594  MVT::MvAddMv( one, *newV, zero, *newV, *lclV );
595 
596  // done with local pointers
597  lclV = Teuchos::null;
598  }
599 
600  // put data into z_, make sure old information is not still hanging around.
601  if (newstate.z != z_) {
602  z_->putScalar();
603  Teuchos::SerialDenseMatrix<int,ScalarType> newZ(Teuchos::View,*newstate.z,curDim_+blockSize_,blockSize_);
605  lclZ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*z_,curDim_+blockSize_,blockSize_) );
606  lclZ->assign(newZ);
607 
608  // done with local pointers
609  lclZ = Teuchos::null;
610  }
611 
612  }
613  else {
614 
615  TEUCHOS_TEST_FOR_EXCEPTION(newstate.V == Teuchos::null,std::invalid_argument,
616  "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial kernel V_0.");
617 
618  TEUCHOS_TEST_FOR_EXCEPTION(newstate.z == Teuchos::null,std::invalid_argument,
619  "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial norms z_0.");
620  }
621 
622  // the solver is initialized
623  initialized_ = true;
624 
625  if (om_->isVerbosity( Debug ) ) {
626  // Check almost everything here
627  CheckList chk;
628  chk.checkV = true;
629  chk.checkArn = true;
630  om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
631  }
632 
633  }
634 
635 
637  // Iterate until the status test informs us we should stop.
638  template <class ScalarType, class MV, class OP>
640  {
641  //
642  // Allocate/initialize data structures
643  //
644  if (initialized_ == false) {
645  initialize();
646  }
647 
648  // Compute the current search dimension.
649  int searchDim = blockSize_*numBlocks_;
650 
652  // iterate until the status test tells us to stop.
653  //
654  // also break if our basis is full
655  //
656  while (stest_->checkStatus(this) != Passed && curDim_+blockSize_ <= searchDim) {
657 
658  iter_++;
659 
660  // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
661  int lclDim = curDim_ + blockSize_;
662 
663  // Get the current part of the basis.
664  std::vector<int> curind(blockSize_);
665  for (int i=0; i<blockSize_; i++) { curind[i] = lclDim + i; }
666  Teuchos::RCP<MV> Vnext = MVT::CloneViewNonConst(*V_,curind);
667 
668  // Get a view of the previous vectors.
669  // This is used for orthogonalization and for computing V^H K H
670  for (int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
671  Teuchos::RCP<const MV> Vprev = MVT::CloneView(*V_,curind);
672 
673  // Compute the next std::vector in the Krylov basis: Vnext = Op*Vprev
674  lp_->apply(*Vprev,*Vnext);
675  Vprev = Teuchos::null;
676 
677  // Remove all previous Krylov basis vectors from Vnext
678  // Get a view of all the previous vectors
679  std::vector<int> prevind(lclDim);
680  for (int i=0; i<lclDim; i++) { prevind[i] = i; }
681  Vprev = MVT::CloneView(*V_,prevind);
682  Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
683 
684  // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
687  ( Teuchos::View,*H_,lclDim,blockSize_,0,curDim_ ) );
689  AsubH.append( subH );
690 
691  // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
694  ( Teuchos::View,*H_,blockSize_,blockSize_,lclDim,curDim_ ) );
695  subH2->putScalar(); // Initialize subdiagonal to zero
696  int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subH2,AVprev);
697 
698  // Copy over the coefficients if we are saving the upper Hessenberg matrix,
699  // just in case we run into an error.
700  if (keepHessenberg_) {
701  // Copy over the orthogonalization coefficients.
704  ( Teuchos::View,*R_,lclDim,blockSize_,0,curDim_ ) );
705  subR->assign(*subH);
706 
707  // Copy over the lower diagonal block of the Hessenberg matrix.
710  ( Teuchos::View,*R_,blockSize_,blockSize_,lclDim,curDim_ ) );
711  subR2->assign(*subH2);
712  }
713 
715  "Belos::BlockGmresIter::iterate(): couldn't generate basis of full rank.");
716  //
717  // V has been extended, and H has been extended.
718  //
719  // Update the QR factorization of the upper Hessenberg matrix
720  //
721  updateLSQR();
722  //
723  // Update basis dim and release all pointers.
724  //
725  Vnext = Teuchos::null;
726  curDim_ += blockSize_;
727  //
728  // When required, monitor some orthogonalities
729  if (om_->isVerbosity( Debug ) ) {
730  // Check almost everything here
731  CheckList chk;
732  chk.checkV = true;
733  chk.checkArn = true;
734  om_->print( Debug, accuracyCheck(chk, ": after local update") );
735  }
736  else if (om_->isVerbosity( OrthoDetails ) ) {
737  CheckList chk;
738  chk.checkV = true;
739  om_->print( OrthoDetails, accuracyCheck(chk, ": after local update") );
740  }
741 
742  } // end while (statusTest == false)
743 
744  }
745 
746 
747  template<class ScalarType, class MV, class OP>
749  {
750  int i, j, maxidx;
751  ScalarType sigma, mu, vscale, maxelem;
752  const ScalarType zero = SCT::zero();
753 
754  // Get correct dimension based on input "dim"
755  // Remember that ortho failures result in an exit before updateLSQR() is called.
756  // Therefore, it is possible that dim == curDim_.
757  int curDim = curDim_;
758  if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
759  curDim = dim;
760  }
761 
763  //
764  // Apply previous transformations and compute new transformation to reduce upper-Hessenberg
765  // system to upper-triangular form.
766  //
767  if (blockSize_ == 1) {
768  //
769  // QR factorization of Least-Squares system with Givens rotations
770  //
771  for (i=0; i<curDim; i++) {
772  //
773  // Apply previous Givens rotations to new column of Hessenberg matrix
774  //
775  blas.ROT( 1, &(*R_)(i,curDim), 1, &(*R_)(i+1, curDim), 1, &cs[i], &sn[i] );
776  }
777  //
778  // Calculate new Givens rotation
779  //
780  blas.ROTG( &(*R_)(curDim,curDim), &(*R_)(curDim+1,curDim), &cs[curDim], &sn[curDim] );
781  (*R_)(curDim+1,curDim) = zero;
782  //
783  // Update RHS w/ new transformation
784  //
785  blas.ROT( 1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim] );
786  }
787  else {
788  //
789  // QR factorization of Least-Squares system with Householder reflectors
790  //
791  for (j=0; j<blockSize_; j++) {
792  //
793  // Apply previous Householder reflectors to new block of Hessenberg matrix
794  //
795  for (i=0; i<curDim+j; i++) {
796  sigma = blas.DOT( blockSize_, &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
797  sigma += (*R_)(i,curDim+j);
798  sigma *= SCT::conjugate(beta[i]);
799  blas.AXPY(blockSize_, ScalarType(-sigma), &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
800  (*R_)(i,curDim+j) -= sigma;
801  }
802  //
803  // Compute new Householder reflector
804  //
805  maxidx = blas.IAMAX( blockSize_+1, &(*R_)(curDim+j,curDim+j), 1 );
806  maxelem = SCT::magnitude((*R_)(curDim+j+maxidx-1,curDim+j));
807  for (i=0; i<blockSize_+1; i++)
808  (*R_)(curDim+j+i,curDim+j) /= maxelem;
809  sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j), 1,
810  &(*R_)(curDim+j+1,curDim+j), 1 );
811  MagnitudeType sign_Rjj = -SCT::real((*R_)(curDim+j,curDim+j)) /
812  SCT::magnitude(SCT::real(((*R_)(curDim+j,curDim+j))));
813  if (sigma == zero) {
814  beta[curDim + j] = zero;
815  } else {
816  mu = SCT::squareroot(SCT::conjugate((*R_)(curDim+j,curDim+j))*(*R_)(curDim+j,curDim+j)+sigma);
817  vscale = (*R_)(curDim+j,curDim+j) - Teuchos::as<ScalarType>(sign_Rjj)*mu;
818  beta[curDim+j] = -Teuchos::as<ScalarType>(sign_Rjj) * vscale / mu;
819  (*R_)(curDim+j,curDim+j) = Teuchos::as<ScalarType>(sign_Rjj)*maxelem*mu;
820  for (i=0; i<blockSize_; i++)
821  (*R_)(curDim+j+1+i,curDim+j) /= vscale;
822  }
823  //
824  // Apply new Householder reflector to rhs
825  //
826  for (i=0; i<blockSize_; i++) {
827  sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j),
828  1, &(*z_)(curDim+j+1,i), 1);
829  sigma += (*z_)(curDim+j,i);
830  sigma *= SCT::conjugate(beta[curDim+j]);
831  blas.AXPY(blockSize_, ScalarType(-sigma), &(*R_)(curDim+j+1,curDim+j),
832  1, &(*z_)(curDim+j+1,i), 1);
833  (*z_)(curDim+j,i) -= sigma;
834  }
835  }
836  } // end if (blockSize_ == 1)
837 
838  // If the least-squares problem is updated wrt "dim" then update the curDim_.
839  if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
840  curDim_ = dim + blockSize_;
841  }
842 
843  } // end updateLSQR()
844 
846  // Check accuracy, orthogonality, and other debugging stuff
847  //
848  // bools specify which tests we want to run (instead of running more than we actually care about)
849  //
850  // checkV : V orthonormal
851  //
852  // checkArn: check the Arnoldi factorization
853  //
854  // NOTE: This method needs to check the current dimension of the subspace, since it is possible to
855  // call this method when curDim_ = 0 (after initialization).
856  template <class ScalarType, class MV, class OP>
857  std::string BlockGmresIter<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const
858  {
859  std::stringstream os;
860  os.precision(2);
861  os.setf(std::ios::scientific, std::ios::floatfield);
862  MagnitudeType tmp;
863 
864  os << " Debugging checks: iteration " << iter_ << where << std::endl;
865 
866  // index vectors for V and F
867  std::vector<int> lclind(curDim_);
868  for (int i=0; i<curDim_; i++) lclind[i] = i;
869  std::vector<int> bsind(blockSize_);
870  for (int i=0; i<blockSize_; i++) { bsind[i] = curDim_ + i; }
871 
872  Teuchos::RCP<const MV> lclV,lclF;
873  Teuchos::RCP<MV> lclAV;
874  if (curDim_)
875  lclV = MVT::CloneView(*V_,lclind);
876  lclF = MVT::CloneView(*V_,bsind);
877 
878  if (chk.checkV) {
879  if (curDim_) {
880  tmp = ortho_->orthonormError(*lclV);
881  os << " >> Error in V^H M V == I : " << tmp << std::endl;
882  }
883  tmp = ortho_->orthonormError(*lclF);
884  os << " >> Error in F^H M F == I : " << tmp << std::endl;
885  if (curDim_) {
886  tmp = ortho_->orthogError(*lclV,*lclF);
887  os << " >> Error in V^H M F == 0 : " << tmp << std::endl;
888  }
889  }
890 
891  if (chk.checkArn) {
892 
893  if (curDim_) {
894  // Compute AV
895  lclAV = MVT::Clone(*V_,curDim_);
896  lp_->apply(*lclV,*lclAV);
897 
898  // Compute AV - VH
899  const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
901  MVT::MvTimesMatAddMv( -one, *lclV, subH, one, *lclAV );
902 
903  // Compute FB_k^T - (AV-VH)
905  blockSize_,curDim_, curDim_ );
906  MVT::MvTimesMatAddMv( -one, *lclF, curB, one, *lclAV );
907 
908  // Compute || FE_k^T - (AV-VH) ||
909  std::vector<MagnitudeType> arnNorms( curDim_ );
910  ortho_->norm( *lclAV, arnNorms );
911 
912  for (int i=0; i<curDim_; i++) {
913  os << " >> Error in Krylov factorization (R = AV-VH-FB^H), ||R[" << i << "]|| : " << arnNorms[i] << std::endl;
914  }
915  }
916  }
917 
918  os << std::endl;
919 
920  return os.str();
921  }
922 
923 } // end Belos namespace
924 
925 #endif /* BELOS_BLOCK_GMRES_ITER_HPP */
ScalarType * values() const
Collection of types and exceptions used within the Belos solvers.
void TRSM(ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
void AXPY(const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
Array< T > & append(const T &x)
Class which manages the output and verbosity of the Belos solvers.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
This class implements the block GMRES iteration, where a block Krylov subspace is constructed...
Teuchos::RCP< const MV > V
The current Krylov basis.
T & get(const std::string &name, T def_value)
MultiVecTraits< ScalarType, MV > MVT
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Declaration of basic traits for the multivector type.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
OrdinalType IAMAX(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Structure to contain pointers to GmresIteration state variables.
void resetNumIters(int iter=0)
Reset the iteration count.
SCT::magnitudeType MagnitudeType
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear prob...
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
virtual ~BlockGmresIter()
Destructor.
Traits class which defines basic operations on multivectors.
bool isParameter(const std::string &name) const
bool isInitialized()
States whether the solver has been initialized or not.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int curDim
The current dimension of the reduction.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
ScalarType DOT(const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
Teuchos::ScalarTraits< ScalarType > SCT
GmresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
OrdinalType numCols() const
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
Class which defines basic traits for the operator type.
void initializeGmres(GmresIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
BlockGmresIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
BlockGmresIter constructor with linear problem, solver utilities, and parameter list of solver option...
Belos&#39;s templated virtual class for providing routines for orthogonalization and orthonormzalition of...
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem...
void setBlockSize(int blockSize)
Set the blocksize.
OrdinalType stride() const
OrdinalType numRows() const
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
int getNumIters() const
Get the current iteration count.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.

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