42 #ifndef BELOS_CG_SINGLE_RED_ITER_HPP
43 #define BELOS_CG_SINGLE_RED_ITER_HPP
77 template<
class ScalarType,
class MV,
class OP>
200 "Belos::CGSingleRedIter::setBlockSize(): Cannot use a block size that is not one.");
249 bool stateStorageInitialized_;
280 template<
class ScalarType,
class MV,
class OP>
289 stateStorageInitialized_(false),
296 template <
class ScalarType,
class MV,
class OP>
299 if (!stateStorageInitialized_) {
304 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
305 stateStorageInitialized_ =
false;
312 if (R_ == Teuchos::null) {
316 "Belos::CGSingleRedIter::setStateSize(): linear problem does not specify multivectors to clone from.");
317 S_ = MVT::Clone( *tmp, 2 );
318 Z_ = MVT::Clone( *tmp, 1 );
319 P_ = MVT::Clone( *tmp, 1 );
320 AP_ = MVT::Clone( *tmp, 1 );
323 std::vector<int> index(1,0);
324 R_ = MVT::CloneViewNonConst( *S_, index );
326 AZ_ = MVT::CloneViewNonConst( *S_, index );
330 stateStorageInitialized_ =
true;
338 template <
class ScalarType,
class MV,
class OP>
342 if (!stateStorageInitialized_)
346 "Belos::CGSingleRedIter::initialize(): Cannot initialize state storage!");
350 std::string errstr(
"Belos::CGSingleRedIter::initialize(): Specified multivectors must have a consistent length and width.");
352 if (newstate.
R != Teuchos::null) {
355 std::invalid_argument, errstr );
357 std::invalid_argument, errstr );
360 if (newstate.
R != R_) {
362 MVT::Assign( *newstate.
R, *R_ );
368 if ( lp_->getLeftPrec() != Teuchos::null ) {
369 lp_->applyLeftPrec( *R_, *Z_ );
370 if ( lp_->getRightPrec() != Teuchos::null ) {
372 lp_->applyRightPrec( *Z_, *tmp );
376 else if ( lp_->getRightPrec() != Teuchos::null ) {
377 lp_->applyRightPrec( *R_, *Z_ );
382 MVT::Assign( *Z_, *P_ );
385 lp_->applyOp( *Z_, *AZ_ );
388 MVT::Assign( *AZ_, *AP_ );
393 "Belos::CGSingleRedIter::initialize(): CGIterationState does not have initial residual.");
403 template <
class ScalarType,
class MV,
class OP>
409 if (initialized_ ==
false) {
415 ScalarType rHz, rHz_old, alpha, beta, delta;
426 "Belos::CGSingleRedIter::iterate(): current linear system has more than one vector!" );
429 MVT::MvTransMv( one, *S_, *Z_, sHz );
436 "Belos::CGSingleRedIter::iterate(): non-positive value for p^H*A*p encountered!" );
445 MVT::MvAddMv( one, *cur_soln_vec, alpha, *P_, *cur_soln_vec );
446 lp_->updateSolution();
450 MVT::MvAddMv( one, *R_, -alpha, *AP_, *R_ );
454 if (stest_->checkStatus(
this) ==
Passed) {
462 if ( lp_->getLeftPrec() != Teuchos::null ) {
463 lp_->applyLeftPrec( *R_, *Z_ );
464 if ( lp_->getRightPrec() != Teuchos::null ) {
466 lp_->applyRightPrec( *Z_, *tmp );
470 else if ( lp_->getRightPrec() != Teuchos::null ) {
471 lp_->applyRightPrec( *R_, *Z_ );
478 lp_->applyOp( *Z_, *AZ_ );
481 MVT::MvTransMv( one, *S_, *Z_, sHz );
488 beta = rHz / rHz_old;
489 alpha = rHz / (delta - (beta*rHz / alpha));
493 "Belos::CGSingleRedIter::iterate(): non-positive value for p^H*A*p encountered!" );
497 MVT::MvAddMv( one, *Z_, beta, *P_, *P_ );
503 MVT::MvAddMv( one, *AZ_, beta, *AP_, *AP_ );