10#ifndef BELOS_PSEUDO_BLOCK_GMRES_ITER_HPP
11#define BELOS_PSEUDO_BLOCK_GMRES_ITER_HPP
28#include "Teuchos_BLAS.hpp"
29#include "Teuchos_SerialDenseMatrix.hpp"
30#include "Teuchos_SerialDenseVector.hpp"
31#include "Teuchos_ScalarTraits.hpp"
32#include "Teuchos_ParameterList.hpp"
33#include "Teuchos_TimeMonitor.hpp"
57 template <
class ScalarType,
class MV>
60 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
69 std::vector<Teuchos::RCP<const MV> >
V;
75 std::vector<Teuchos::RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > >
H;
77 std::vector<Teuchos::RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > >
R;
79 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,ScalarType> > >
Z;
81 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,ScalarType> > >
sn;
82 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,MagnitudeType> > >
cs;
106 template<
class ScalarType,
class MV,
class OP>
116 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
134 Teuchos::ParameterList ¶ms );
209 state.
V.resize(numRHS_);
210 state.
H.resize(numRHS_);
211 state.
Z.resize(numRHS_);
212 state.
sn.resize(numRHS_);
213 state.
cs.resize(numRHS_);
214 for (
int i=0; i<numRHS_; ++i) {
218 state.
sn[i] = sn_[i];
219 state.
cs[i] = cs_[i];
271 if (!initialized_)
return 0;
292 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
293 "Belos::PseudoBlockGmresIter::setBlockSize(): Cannot use a block size that is not one.");
312 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
313 const Teuchos::RCP<OutputManager<ScalarType> > om_;
314 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
315 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
326 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > > sn_;
327 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,MagnitudeType> > > cs_;
330 Teuchos::RCP<MV> U_vec_, AU_vec_;
333 Teuchos::RCP<MV> cur_block_rhs_, cur_block_sol_;
349 std::vector<Teuchos::RCP<MV> > V_;
354 std::vector<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > H_;
359 std::vector<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > R_;
360 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > > Z_;
365 template<
class ScalarType,
class MV,
class OP>
370 Teuchos::ParameterList ¶ms ):
382 TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter(
"Num Blocks"), std::invalid_argument,
383 "Belos::PseudoBlockGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
384 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
391 template <
class ScalarType,
class MV,
class OP>
397 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0, std::invalid_argument,
"Belos::PseudoBlockGmresIter::setNumBlocks was passed a non-positive argument.");
399 numBlocks_ = numBlocks;
402 initialized_ =
false;
407 template <
class ScalarType,
class MV,
class OP>
414 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
416 return currentUpdate;
418 currentUpdate =
MVT::Clone(*(V_[0]), numRHS_);
419 std::vector<int> index(1), index2(curDim_);
420 for (
int i=0; i<curDim_; ++i) {
423 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
424 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
425 Teuchos::BLAS<int,ScalarType> blas;
427 for (
int i=0; i<numRHS_; ++i) {
433 Teuchos::SerialDenseVector<int,ScalarType> y( Teuchos::Copy, Z_[i]->values(), curDim_ );
437 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
438 Teuchos::NON_UNIT_DIAG, curDim_, 1, one,
439 H_[i]->values(), H_[i]->stride(), y.values(), y.stride() );
445 return currentUpdate;
452 template <
class ScalarType,
class MV,
class OP>
453 Teuchos::RCP<const MV>
457 typedef typename Teuchos::ScalarTraits<ScalarType> STS;
463 if (
static_cast<int> (norms->size()) < numRHS_)
464 norms->resize (numRHS_);
466 Teuchos::BLAS<int, ScalarType> blas;
467 for (
int j = 0; j < numRHS_; ++j)
469 const ScalarType curNativeResid = (*Z_[j])(curDim_);
470 (*norms)[j] = STS::magnitude (curNativeResid);
473 return Teuchos::null;
477 template <
class ScalarType,
class MV,
class OP>
492 std::string errstr (
"Belos::PseudoBlockGmresIter::initialize(): "
493 "Specified multivectors must have a consistent "
494 "length and width.");
497 TEUCHOS_TEST_FOR_EXCEPTION((
int)newstate.
V.size()==0 || (
int)newstate.
Z.size()==0,
498 std::invalid_argument,
499 "Belos::PseudoBlockGmresIter::initialize(): "
500 "V and/or Z was not specified in the input state; "
501 "the V and/or Z arrays have length zero.");
512 RCP<const MV> lhsMV = lp_->getLHS();
513 RCP<const MV> rhsMV = lp_->getRHS();
517 RCP<const MV> vectorInBasisSpace = rhsMV.is_null() ? lhsMV : rhsMV;
520 TEUCHOS_TEST_FOR_EXCEPTION(vectorInBasisSpace.is_null(),
521 std::invalid_argument,
522 "Belos::PseudoBlockGmresIter::initialize(): "
523 "The linear problem to solve does not specify multi"
524 "vectors from which we can clone basis vectors. The "
525 "right-hand side(s), left-hand side(s), or both should "
530 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
curDim > numBlocks_+1,
531 std::invalid_argument,
533 curDim_ = newstate.
curDim;
539 for (
int i=0; i<numRHS_; ++i) {
544 V_[i] =
MVT::Clone (*vectorInBasisSpace, numBlocks_ + 1);
549 std::invalid_argument, errstr );
551 std::invalid_argument, errstr );
557 if (newstate.
V[i] != V_[i]) {
559 if (curDim_ == 0 && lclDim > 1) {
561 <<
"Belos::PseudoBlockGmresIter::initialize(): the solver was "
562 <<
"initialized with a kernel of " << lclDim
564 <<
"The block size however is only " << 1
566 <<
"The last " << lclDim - 1 <<
" vectors will be discarded."
569 std::vector<int> nevind (curDim_ + 1);
570 for (
int j = 0; j < curDim_ + 1; ++j)
575 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
576 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
580 lclV = Teuchos::null;
587 for (
int i=0; i<numRHS_; ++i) {
589 if (Z_[i] == Teuchos::null) {
590 Z_[i] = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,ScalarType>() );
592 if (Z_[i]->length() < numBlocks_+1) {
593 Z_[i]->shapeUninitialized(numBlocks_+1, 1);
597 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
Z[i]->numRows() < curDim_, std::invalid_argument, errstr);
600 if (newstate.
Z[i] != Z_[i]) {
604 Teuchos::SerialDenseVector<int,ScalarType> newZ(Teuchos::View,newstate.
Z[i]->values(),curDim_+1);
605 Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > lclZ;
606 lclZ = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,ScalarType>(Teuchos::View,Z_[i]->values(),curDim_+1) );
610 lclZ = Teuchos::null;
617 for (
int i=0; i<numRHS_; ++i) {
619 if (H_[i] == Teuchos::null) {
620 H_[i] = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>() );
622 if (H_[i]->numRows() < numBlocks_+1 || H_[i]->numCols() < numBlocks_) {
623 H_[i]->shapeUninitialized(numBlocks_+1, numBlocks_);
627 if ((
int)newstate.
H.size() == numRHS_) {
630 TEUCHOS_TEST_FOR_EXCEPTION((newstate.
H[i]->numRows() < curDim_ || newstate.
H[i]->numCols() < curDim_), std::invalid_argument,
631 "Belos::PseudoBlockGmresIter::initialize(): Specified Hessenberg matrices must have a consistent size to the current subspace dimension");
633 if (newstate.
H[i] != H_[i]) {
636 Teuchos::SerialDenseMatrix<int,ScalarType> newH(Teuchos::View,*newstate.
H[i],curDim_+1, curDim_);
637 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclH;
638 lclH = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*H_[i],curDim_+1, curDim_) );
642 lclH = Teuchos::null;
654 if ((
int)newstate.
cs.size() == numRHS_ && (
int)newstate.
sn.size() == numRHS_) {
655 for (
int i=0; i<numRHS_; ++i) {
656 if (cs_[i] != newstate.
cs[i])
657 cs_[i] = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,MagnitudeType>(*newstate.
cs[i]) );
658 if (sn_[i] != newstate.
sn[i])
659 sn_[i] = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,ScalarType>(*newstate.
sn[i]) );
664 for (
int i=0; i<numRHS_; ++i) {
665 if (cs_[i] == Teuchos::null)
666 cs_[i] = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,MagnitudeType>(numBlocks_+1) );
668 cs_[i]->resize(numBlocks_+1);
669 if (sn_[i] == Teuchos::null)
670 sn_[i] = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,ScalarType>(numBlocks_+1) );
672 sn_[i]->resize(numBlocks_+1);
683 template <
class ScalarType,
class MV,
class OP>
689 if (initialized_ ==
false) {
693 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
694 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
697 int searchDim = numBlocks_;
702 std::vector<int> index(1);
703 std::vector<int> index2(1);
705 Teuchos::RCP<MV> U_vec =
MVT::Clone( *V_[0], numRHS_ );
708 Teuchos::RCP<MV> AU_vec =
MVT::Clone( *V_[0], numRHS_ );
710 for (
int i=0; i<numRHS_; ++i) {
714 MVT::MvAddMv( one, *tmp_vec, zero, *tmp_vec, *U_vec_view );
722 while (stest_->checkStatus(
this) !=
Passed && curDim_ < searchDim) {
728 lp_->apply( *U_vec, *AU_vec );
733 int num_prev = curDim_+1;
734 index.resize( num_prev );
735 for (
int i=0; i<num_prev; ++i) {
741 for (
int i=0; i<numRHS_; ++i) {
746 Teuchos::Array< Teuchos::RCP<const MV> > V_array( 1, V_prev );
755 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > h_new
756 = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
757 ( Teuchos::View, *H_[i], num_prev, 1, 0, curDim_ ) );
758 Teuchos::Array< Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > h_array( 1, h_new );
760 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > r_new
761 = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
762 ( Teuchos::View, *H_[i], 1, 1, num_prev, curDim_ ) );
767 ortho_->projectAndNormalize( *V_new, h_array, r_new, V_array );
772 index2[0] = curDim_+1;
780 Teuchos::RCP<MV> tmp_AU_vec = U_vec;
800 template<
class ScalarType,
class MV,
class OP>
806 int curDim = curDim_;
812 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
814 Teuchos::BLAS<int, ScalarType> blas;
816 for (i=0; i<numRHS_; ++i) {
822 for (j=0; j<curDim; j++) {
826 blas.ROT( 1, &(*H_[i])(j,curDim), 1, &(*H_[i])(j+1, curDim), 1, &(*cs_[i])[j], &(*sn_[i])[j] );
831 blas.ROTG( &(*H_[i])(curDim,curDim), &(*H_[i])(curDim+1,curDim), &(*cs_[i])[curDim], &(*sn_[i])[curDim] );
832 (*H_[i])(curDim+1,curDim) = zero;
836 blas.ROT( 1, &(*Z_[i])(curDim), 1, &(*Z_[i])(curDim+1), 1, &(*cs_[i])[curDim], &(*sn_[i])[curDim] );
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which describes the basic interface to the linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
BelosError(const std::string &what_arg)
Iteration()
Default Constructor.
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
Traits class which defines basic operations on multivectors.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
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).
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).
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Class which defines basic traits for the operator type.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
PseudoBlockGmresIterOrthoFailure(const std::string &what_arg)
MultiVecTraits< ScalarType, MV > MVT
bool isInitialized()
States whether the solver has been initialized or not.
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
SCT::magnitudeType MagnitudeType
OperatorTraits< ScalarType, MV, OP > OPT
PseudoBlockGmresIter(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 ¶ms)
PseudoBlockGmresIter constructor with linear problem, solver utilities, and parameter list of solver ...
void initialize(const PseudoBlockGmresIterState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem.
Teuchos::ScalarTraits< ScalarType > SCT
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the "native" residual vectors.
void setBlockSize(int blockSize)
Set the blocksize.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
void resetNumIters(int iter=0)
Reset the iteration count.
PseudoBlockGmresIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
int getNumIters() const
Get the current iteration count.
virtual ~PseudoBlockGmresIter()
Destructor.
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to PseudoBlockGmresIter state variables.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > sn
The current Given's rotation coefficients.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > R
The current upper-triangular matrix from the QR reduction of H.
std::vector< Teuchos::RCP< const MV > > V
The current Krylov basis.
PseudoBlockGmresIterState()
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > Z
The current right-hand side of the least squares system RY = Z.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > H
The current Hessenberg matrix.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, MagnitudeType > > > cs
Teuchos::ScalarTraits< ScalarType > SCT
int curDim
The current dimension of the reduction.
SCT::magnitudeType MagnitudeType