42 #ifndef BELOS_PSEUDO_BLOCK_CG_ITER_MP_VECTOR_HPP 43 #define BELOS_PSEUDO_BLOCK_CG_ITER_MP_VECTOR_HPP 54 #include "BelosPseudoBlockCGIter.hpp" 73 template<
class Storage,
class MV,
class OP>
75 virtual public CGIteration<Sacado::MP::Vector<Storage>, MV, OP> {
83 typedef MultiVecTraits<ScalarType,MV>
MVT;
84 typedef OperatorTraits<ScalarType,MV,OP>
OPT;
99 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
145 void initializeCG(CGIterationState<ScalarType,MV>& newstate);
152 CGIterationState<ScalarType,MV> empty;
164 CGIterationState<ScalarType,MV> state;
199 const LinearProblem<ScalarType,MV,OP>&
getProblem()
const {
return *lp_; }
207 "Belos::PseudoBlockCGIter::setBlockSize(): Cannot use a block size that is not one.");
224 if (static_cast<size_type> (iter_) >= diag_.size ()) {
227 return diag_ (0, iter_);
239 if (static_cast<size_type> (iter_) >= offdiag_.size ()) {
242 return offdiag_ (0, iter_);
302 template<
class StorageType,
class MV,
class OP>
304 PseudoBlockCGIter(
const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
305 const Teuchos::RCP<OutputManager<ScalarType> > &printer,
306 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
314 assertPositiveDefiniteness_( params.get(
"Assert Positive Definiteness",
true) ),
315 numEntriesForCondEst_(params.get(
"Max Size For Condest",0) ),
316 breakDownTol_(params.get(
"Ensemble Breakdown Tolerance", 0.0))
323 template <
class StorageType,
class MV,
class OP>
324 void PseudoBlockCGIter<Sacado::MP::Vector<StorageType>,MV,OP>::
325 initializeCG(CGIterationState<ScalarType,MV>& newstate)
331 "Belos::PseudoBlockCGIter::initialize(): Cannot initialize state storage!");
337 int numRHS = MVT::GetNumberVecs(*tmp);
343 R_ = MVT::Clone( *tmp, numRHS_ );
344 Z_ = MVT::Clone( *tmp, numRHS_ );
345 P_ = MVT::Clone( *tmp, numRHS_ );
346 AP_ = MVT::Clone( *tmp, numRHS_ );
350 if(numEntriesForCondEst_ > 0) {
351 diag_.resize(numEntriesForCondEst_);
352 offdiag_.resize(numEntriesForCondEst_-1);
357 std::string errstr(
"Belos::BlockPseudoCGIter::initialize(): Specified multivectors must have a consistent length and width.");
366 std::invalid_argument, errstr );
368 std::invalid_argument, errstr );
371 if (newstate.R != R_) {
373 MVT::MvAddMv( one, *newstate.R, zero, *newstate.R, *R_ );
379 if ( lp_->getLeftPrec() != Teuchos::null ) {
380 lp_->applyLeftPrec( *R_, *Z_ );
381 if ( lp_->getRightPrec() != Teuchos::null ) {
383 lp_->applyRightPrec( *Z_, *tmp1 );
387 else if ( lp_->getRightPrec() != Teuchos::null ) {
388 lp_->applyRightPrec( *R_, *Z_ );
393 MVT::MvAddMv( one, *Z_, zero, *Z_, *P_ );
398 "Belos::CGIter::initialize(): CGStateIterState does not have initial residual.");
408 template <
class StorageType,
class MV,
class OP>
409 void PseudoBlockCGIter<Sacado::MP::Vector<StorageType>,MV,OP>::
415 if (initialized_ ==
false) {
421 std::vector<int> index(1);
422 std::vector<ScalarType> rHz( numRHS_ ), rHz_old( numRHS_ ), pAp( numRHS_ );
430 ScalarType pAp_old=one, beta_old=one ,rHz_old2=one;
436 MVT::MvDot( *R_, *Z_, rHz );
438 if ( assertPositiveDefiniteness_ )
439 for (i=0; i<numRHS_; ++i)
442 "Belos::PseudoBlockCGIter::iterate(): negative value for r^H*M*r encountered!" );
447 while (stest_->checkStatus(
this) != Passed) {
453 lp_->applyOp( *P_, *AP_ );
456 MVT::MvDot( *P_, *AP_, pAp );
458 for (i=0; i<numRHS_; ++i) {
459 if ( assertPositiveDefiniteness_ )
463 "Belos::PseudoBlockCGIter::iterate(): non-positive value for p^H*A*p encountered!" );
465 const int ensemble_size = pAp[i].size();
466 for (
int k=0; k<ensemble_size; ++k) {
467 if (SVT::magnitude(pAp[i].coeff(k)) <= breakDownTol_)
468 alpha(i,i).fastAccessCoeff(k) = 0.0;
470 alpha(i,i).fastAccessCoeff(k) = rHz[i].coeff(k) / pAp[i].coeff(k);
477 MVT::MvTimesMatAddMv( one, *P_, alpha, one, *cur_soln_vec );
478 lp_->updateSolution();
482 for (i=0; i<numRHS_; ++i) {
488 MVT::MvTimesMatAddMv( -one, *AP_, alpha, one, *R_ );
493 if ( lp_->getLeftPrec() != Teuchos::null ) {
494 lp_->applyLeftPrec( *R_, *Z_ );
495 if ( lp_->getRightPrec() != Teuchos::null ) {
497 lp_->applyRightPrec( *Z_, *tmp );
501 else if ( lp_->getRightPrec() != Teuchos::null ) {
502 lp_->applyRightPrec( *R_, *Z_ );
508 MVT::MvDot( *R_, *Z_, rHz );
509 if ( assertPositiveDefiniteness_ )
510 for (i=0; i<numRHS_; ++i)
513 "Belos::PseudoBlockCGIter::iterate(): negative value for r^H*M*r encountered!" );
516 for (i=0; i<numRHS_; ++i) {
517 const int ensemble_size = rHz[i].size();
518 for (
int k=0; k<ensemble_size; ++k) {
519 if (SVT::magnitude(rHz[i].coeff(k)) <= breakDownTol_)
520 beta(i,i).fastAccessCoeff(k) = 0.0;
522 beta(i,i).fastAccessCoeff(k) = rHz[i].coeff(k) / rHz_old[i].coeff(k);
527 MVT::MvAddMv( one, *Z_i, beta(i,i), *P_i, *P_i );
539 rHz_old2 = rHz_old[0];
540 beta_old = beta(0,0);
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
bool assertPositiveDefiniteness_
Teuchos::ScalarTraits< typename Storage::value_type > SVT
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
static magnitudeType real(T a)
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
OperatorTraits< ScalarType, MV, OP > OPT
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
virtual ~PseudoBlockCGIter()
Destructor.
int getNumIters() const
Get the current iteration count.
This class implements the pseudo-block CG iteration, where the basic CG algorithm is performed on all...
Teuchos::ScalarTraits< ScalarType > SCT
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void setBlockSize(int blockSize)
Set the blocksize.
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation.
MultiVecTraits< ScalarType, MV > MVT
int numEntriesForCondEst_
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation.
SCT::magnitudeType MagnitudeType
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
Sacado::MP::Vector< Storage > ScalarType
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
bool isInitialized()
States whether the solver has been initialized or not.
SVT::magnitudeType breakDownTol_
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
const Teuchos::RCP< OutputManager< ScalarType > > om_
Teuchos::ArrayRCP< MagnitudeType > offdiag_
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_