44 #ifndef ANASAZI_IRTR_HPP
45 #define ANASAZI_IRTR_HPP
85 template <
class ScalarType,
class MV,
class OP>
151 std::vector<std::string> stopReasons_;
155 const MagnitudeType ZERO;
156 const MagnitudeType ONE;
161 void solveTRSubproblem();
164 MagnitudeType rho_prime_;
167 MagnitudeType normgradf0_;
170 trRetType innerStop_;
173 int innerIters_, totalInnerIters_;
184 template <
class ScalarType,
class MV,
class OP>
193 RTRBase<ScalarType,MV,OP>(problem,sorter,printer,tester,ortho,params,
"IRTR",false),
199 stopReasons_.push_back(
"n/a");
200 stopReasons_.push_back(
"maximum iterations");
201 stopReasons_.push_back(
"negative curvature");
202 stopReasons_.push_back(
"exceeded TR");
203 stopReasons_.push_back(
"kappa convergence");
204 stopReasons_.push_back(
"theta convergence");
206 rho_prime_ = params.
get(
"Rho Prime",0.5);
208 "Anasazi::IRTR::constructor: rho_prime must be in (0,1).");
210 useSA_ = params.
get<
bool>(
"Use SA",
false);
223 template <
class ScalarType,
class MV,
class OP>
234 using Teuchos::tuple;
236 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
243 innerStop_ = MAXIMUM_ITERATIONS;
245 const int n = MVT::GetGlobalLength(*this->eta_);
246 const int p = this->blockSize_;
247 const int d = n*p - (p*p+p)/2;
261 const double D2 = ONE/rho_prime_ - ONE;
263 std::vector<MagnitudeType> d_Hd(p), alpha(p), beta(p), z_r(p), zold_rold(p);
264 std::vector<MagnitudeType> eBe(p), eBd(p), dBd(p), new_eBe(p);
265 MagnitudeType r0_norm;
267 MVT::MvInit(*this->eta_ ,0.0);
268 MVT::MvInit(*this->Aeta_,0.0);
270 MVT::MvInit(*this->Beta_,0.0);
286 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
287 TimeMonitor lcltimer( *this->timerOrtho_ );
289 this->orthman_->projectGen(
291 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
293 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
294 if (this->leftMost_) {
295 MVT::MvScale(*this->R_,2.0);
298 MVT::MvScale(*this->R_,-2.0);
302 r0_norm = MAT::squareroot( RTRBase<ScalarType,MV,OP>::ginner(*this->R_) );
307 MagnitudeType kconv = r0_norm * this->conv_kappa_;
310 MagnitudeType tconv = MAT::pow(r0_norm,this->conv_theta_+ONE);
311 if (this->om_->isVerbosity(
Debug)) {
312 this->om_->stream(
Debug)
313 <<
" >> |r0| : " << r0_norm << endl
314 <<
" >> kappa conv : " << kconv << endl
315 <<
" >> theta conv : " << tconv << endl;
323 if (this->hasPrec_ && this->olsenPrec_)
325 std::vector<int> ind(this->blockSize_);
326 for (
int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
329 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
330 TimeMonitor prectimer( *this->timerPrec_ );
332 OPT::Apply(*this->Prec_,*this->BX_,*PBX);
333 this->counterPrec_ += this->blockSize_;
344 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
345 TimeMonitor prectimer( *this->timerPrec_ );
347 OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
348 this->counterPrec_ += this->blockSize_;
350 if (this->olsenPrec_) {
351 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
352 TimeMonitor orthtimer( *this->timerOrtho_ );
354 this->orthman_->projectGen(
356 tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),
false,
358 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
361 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
362 TimeMonitor orthtimer( *this->timerOrtho_ );
364 this->orthman_->projectGen(
366 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
368 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
370 RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,*this->Z_,z_r);
373 RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,z_r);
376 if (this->om_->isVerbosity(
Debug )) {
378 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
380 if (this->hasPrec_) chk.checkZ =
true;
381 this->om_->print(
Debug, this->accuracyCheck(chk,
"after computing gradient") );
385 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
387 if (this->hasPrec_) chk.checkZ =
true;
388 this->om_->print(
OrthoDetails, this->accuracyCheck(chk,
"after computing gradient") );
392 MVT::MvAddMv(-ONE,*this->Z_,ZERO,*this->Z_,*this->delta_);
394 if (this->om_->isVerbosity(
Debug)) {
395 RCP<MV> Heta = MVT::Clone(*this->eta_,this->blockSize_);
397 MVT::MvAddMv(ONE,*this->Beta_,ZERO,*this->Beta_,*Heta);
398 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
399 MVT::MvScale(*Heta,theta_comp);
400 if (this->leftMost_) {
401 MVT::MvAddMv( 2.0,*this->Aeta_,-2.0,*Heta,*Heta);
404 MVT::MvAddMv(-2.0,*this->Aeta_, 2.0,*Heta,*Heta);
408 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
409 TimeMonitor lcltimer( *this->timerOrtho_ );
411 this->orthman_->projectGen(
413 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
415 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
418 std::vector<MagnitudeType> eg(this->blockSize_),
419 eHe(this->blockSize_);
420 RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*this->AX_,eg);
421 RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*Heta,eHe);
422 if (this->leftMost_) {
423 for (
int j=0; j<this->blockSize_; ++j) {
424 eg[j] = this->theta_[j] + 2*eg[j] + .5*eHe[j];
428 for (
int j=0; j<this->blockSize_; ++j) {
429 eg[j] = -this->theta_[j] - 2*eg[j] + .5*eHe[j];
432 this->om_->stream(
Debug)
433 <<
" Debugging checks: IRTR inner iteration " << innerIters_ << endl
434 <<
" >> m_X(eta) : " << std::accumulate(eg.begin(),eg.end(),0.0) << endl;
435 for (
int j=0; j<this->blockSize_; ++j) {
436 this->om_->stream(
Debug)
437 <<
" >> m_X(eta_" << j <<
") : " << eg[j] << endl;
444 for (innerIters_=1; innerIters_<=d; ++innerIters_) {
452 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
453 TimeMonitor lcltimer( *this->timerAOp_ );
455 OPT::Apply(*this->AOp_,*this->delta_,*this->Adelta_);
456 this->counterAOp_ += this->blockSize_;
459 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
460 TimeMonitor lcltimer( *this->timerBOp_ );
462 OPT::Apply(*this->BOp_,*this->delta_,*this->Bdelta_);
463 this->counterBOp_ += this->blockSize_;
466 MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*this->Bdelta_);
470 RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_ ,*this->Bdelta_,eBd);
471 RTRBase<ScalarType,MV,OP>::ginnersep(*this->delta_,*this->Bdelta_,dBd);
473 MVT::MvAddMv(ONE,*this->Bdelta_,ZERO,*this->Bdelta_,*this->Hdelta_);
475 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
476 MVT::MvScale(*this->Hdelta_,theta_comp);
478 if (this->leftMost_) {
479 MVT::MvAddMv( 2.0,*this->Adelta_,-2.0,*this->Hdelta_,*this->Hdelta_);
482 MVT::MvAddMv(-2.0,*this->Adelta_, 2.0,*this->Hdelta_,*this->Hdelta_);
486 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
487 TimeMonitor lcltimer( *this->timerOrtho_ );
489 this->orthman_->projectGen(
491 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
493 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
495 RTRBase<ScalarType,MV,OP>::ginnersep(*this->delta_,*this->Hdelta_,d_Hd);
499 for (
unsigned int j=0; j<alpha.size(); ++j)
501 alpha[j] = z_r[j]/d_Hd[j];
505 for (
unsigned int j=0; j<alpha.size(); ++j)
507 new_eBe[j] = eBe[j] + 2*alpha[j]*eBd[j] + alpha[j]*alpha[j]*dBd[j];
510 if (this->om_->isVerbosity(
Debug)) {
511 for (
unsigned int j=0; j<alpha.size(); j++) {
512 this->om_->stream(
Debug)
513 <<
" >> z_r[" << j <<
"] : " << z_r[j]
514 <<
" d_Hd[" << j <<
"] : " << d_Hd[j] << endl
515 <<
" >> eBe[" << j <<
"] : " << eBe[j]
516 <<
" neweBe[" << j <<
"] : " << new_eBe[j] << endl
517 <<
" >> eBd[" << j <<
"] : " << eBd[j]
518 <<
" dBd[" << j <<
"] : " << dBd[j] << endl;
523 std::vector<int> trncstep;
527 bool atleastonenegcur =
false;
528 for (
unsigned int j=0; j<d_Hd.size(); ++j) {
530 trncstep.push_back(j);
531 atleastonenegcur =
true;
533 else if (new_eBe[j] > D2) {
534 trncstep.push_back(j);
538 if (!trncstep.empty())
541 if (this->om_->isVerbosity(
Debug)) {
542 for (
unsigned int j=0; j<trncstep.size(); ++j) {
543 this->om_->stream(
Debug)
544 <<
" >> alpha[" << trncstep[j] <<
"] : " << alpha[trncstep[j]] << endl;
547 for (
unsigned int j=0; j<trncstep.size(); ++j) {
548 int jj = trncstep[j];
549 alpha[jj] = ( -eBd[jj] + MAT::squareroot(eBd[jj]*eBd[jj] + dBd[jj]*(D2-eBe[jj]) ) ) / dBd[jj];
551 if (this->om_->isVerbosity(
Debug)) {
552 for (
unsigned int j=0; j<trncstep.size(); ++j) {
553 this->om_->stream(
Debug)
554 <<
" >> tau[" << trncstep[j] <<
"] : " << alpha[trncstep[j]] << endl;
557 if (atleastonenegcur) {
558 innerStop_ = NEGATIVE_CURVATURE;
561 innerStop_ = EXCEEDED_TR;
570 std::vector<ScalarType> alpha_comp(alpha.begin(),alpha.end());
571 MVT::MvScale(*this->delta_,alpha_comp);
572 MVT::MvScale(*this->Adelta_,alpha_comp);
574 MVT::MvScale(*this->Bdelta_,alpha_comp);
576 MVT::MvScale(*this->Hdelta_,alpha_comp);
578 MVT::MvAddMv(ONE,*this->delta_ ,ONE,*this->eta_ ,*this->eta_);
579 MVT::MvAddMv(ONE,*this->Adelta_,ONE,*this->Aeta_,*this->Aeta_);
581 MVT::MvAddMv(ONE,*this->Bdelta_,ONE,*this->Beta_,*this->Beta_);
589 if (this->om_->isVerbosity(
Debug)) {
590 RCP<MV> Heta = MVT::Clone(*this->eta_,this->blockSize_);
592 MVT::MvAddMv(ONE,*this->Beta_,ZERO,*this->Beta_,*Heta);
594 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
595 MVT::MvScale(*Heta,theta_comp);
597 if (this->leftMost_) {
598 MVT::MvAddMv( 2.0,*this->Aeta_,-2.0,*Heta,*Heta);
601 MVT::MvAddMv(-2.0,*this->Aeta_, 2.0,*Heta,*Heta);
605 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
606 TimeMonitor lcltimer( *this->timerOrtho_ );
608 this->orthman_->projectGen(
610 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
612 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
615 std::vector<MagnitudeType> eg(this->blockSize_),
616 eHe(this->blockSize_);
617 RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*this->AX_,eg);
618 RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*Heta,eHe);
619 if (this->leftMost_) {
620 for (
int j=0; j<this->blockSize_; ++j) {
621 eg[j] = this->theta_[j] + 2*eg[j] + .5*eHe[j];
625 for (
int j=0; j<this->blockSize_; ++j) {
626 eg[j] = -this->theta_[j] - 2*eg[j] + .5*eHe[j];
629 this->om_->stream(
Debug)
630 <<
" Debugging checks: IRTR inner iteration " << innerIters_ << endl
631 <<
" >> m_X(eta) : " << std::accumulate(eg.begin(),eg.end(),0.0) << endl;
632 for (
int j=0; j<this->blockSize_; ++j) {
633 this->om_->stream(
Debug)
634 <<
" >> m_X(eta_" << j <<
") : " << eg[j] << endl;
640 if (!trncstep.empty()) {
648 MVT::MvAddMv(ONE,*this->Hdelta_,ONE,*this->R_,*this->R_);
651 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
652 TimeMonitor lcltimer( *this->timerOrtho_ );
654 this->orthman_->projectGen(
656 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
658 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
663 MagnitudeType r_norm = MAT::squareroot(RTRBase<ScalarType,MV,OP>::ginner(*this->R_,*this->R_));
671 if (this->om_->isVerbosity(
Debug)) {
672 this->om_->stream(
Debug)
673 <<
" >> |r" << innerIters_ <<
"| : " << r_norm << endl;
675 if ( r_norm <= ANASAZI_MIN(tconv,kconv) ) {
676 if (tconv <= kconv) {
677 innerStop_ = THETA_CONVERGENCE;
680 innerStop_ = KAPPA_CONVERGENCE;
692 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
693 TimeMonitor prectimer( *this->timerPrec_ );
695 OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
696 this->counterPrec_ += this->blockSize_;
698 if (this->olsenPrec_) {
699 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
700 TimeMonitor orthtimer( *this->timerOrtho_ );
702 this->orthman_->projectGen(
704 tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),
false,
706 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
709 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
710 TimeMonitor orthtimer( *this->timerOrtho_ );
712 this->orthman_->projectGen(
714 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),
false,
716 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
718 RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,*this->Z_,z_r);
721 RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,z_r);
737 for (
unsigned int j=0; j<beta.size(); ++j) {
738 beta[j] = z_r[j]/(zold_rold[j]*alpha[j]);
741 std::vector<ScalarType> beta_comp(beta.begin(),beta.end());
742 MVT::MvScale(*this->delta_,beta_comp);
744 MVT::MvAddMv(-ONE,*this->Z_,ONE,*this->delta_,*this->delta_);
750 if (innerIters_ > d) innerIters_ = d;
752 this->om_->stream(
Debug)
753 <<
" >> stop reason is " << stopReasons_[innerStop_] << endl
761 template <
class ScalarType,
class MV,
class OP>
766 using Teuchos::tuple;
767 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
777 if (this->initialized_ ==
false) {
782 if (useSA_ ==
true) {
783 AA.
reshape(2*this->blockSize_,2*this->blockSize_);
784 BB.
reshape(2*this->blockSize_,2*this->blockSize_);
785 S.
reshape(2*this->blockSize_,2*this->blockSize_);
788 AA.
reshape(this->blockSize_,this->blockSize_);
789 BB.
reshape(this->blockSize_,this->blockSize_);
790 S.
reshape(this->blockSize_,this->blockSize_);
795 innerStop_ = UNINITIALIZED;
798 while (this->tester_->checkStatus(
this) !=
Passed) {
801 if (this->om_->isVerbosity(
Debug)) {
802 this->currentStatus( this->om_->stream(
Debug) );
813 totalInnerIters_ += innerIters_;
816 if (this->om_->isVerbosity(
Debug ) ) {
821 chk.checkAEta =
true;
822 chk.checkBEta =
true;
823 this->om_->print(
Debug, this->accuracyCheck(chk,
"in iterate() after solveTRSubproblem()") );
824 this->om_->stream(
Debug)
829 if (useSA_ ==
false) {
833 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
834 TimeMonitor lcltimer( *this->timerLocalUpdate_ );
836 MVT::MvAddMv(ONE,*this->X_,ONE,*this->eta_,*this->delta_);
837 MVT::MvAddMv(ONE,*this->AX_,ONE,*this->Aeta_,*this->Adelta_);
839 MVT::MvAddMv(ONE,*this->BX_,ONE,*this->Beta_,*this->Bdelta_);
843 if (this->om_->isVerbosity(
Debug ) ) {
848 E(this->blockSize_,this->blockSize_);
849 MVT::MvTransMv(ONE,*this->delta_,*this->Bdelta_,XE);
850 MVT::MvTransMv(ONE,*this->eta_,*this->Beta_,E);
851 this->om_->stream(
Debug)
852 <<
" >> Error in AX+AEta == A(X+Eta) : " << Utils::errorEquality(*this->delta_,*this->Adelta_,this->AOp_) << endl
853 <<
" >> Error in BX+BEta == B(X+Eta) : " << Utils::errorEquality(*this->delta_,*this->Bdelta_,this->BOp_) << endl
854 <<
" >> norm( (X+Eta)^T B (X+Eta) ) : " << XE.normFrobenius() << endl
864 MagnitudeType oldfx = this->fx_;
865 std::vector<MagnitudeType> oldtheta(this->theta_), newtheta(AA.
numRows());
868 if (useSA_ ==
true) {
869 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
870 TimeMonitor lcltimer( *this->timerLocalProj_ );
873 AA12(
Teuchos::View,AA,this->blockSize_,this->blockSize_,0,this->blockSize_),
874 AA22(
Teuchos::View,AA,this->blockSize_,this->blockSize_,this->blockSize_,this->blockSize_);
876 BB12(
Teuchos::View,BB,this->blockSize_,this->blockSize_,0,this->blockSize_),
877 BB22(
Teuchos::View,BB,this->blockSize_,this->blockSize_,this->blockSize_,this->blockSize_);
878 MVT::MvTransMv(ONE,*this->X_ ,*this->AX_ ,AA11);
879 MVT::MvTransMv(ONE,*this->X_ ,*this->Aeta_,AA12);
880 MVT::MvTransMv(ONE,*this->eta_,*this->Aeta_,AA22);
881 MVT::MvTransMv(ONE,*this->X_ ,*this->BX_ ,BB11);
882 MVT::MvTransMv(ONE,*this->X_ ,*this->Beta_,BB12);
883 MVT::MvTransMv(ONE,*this->eta_,*this->Beta_,BB22);
886 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
887 TimeMonitor lcltimer( *this->timerLocalProj_ );
889 MVT::MvTransMv(ONE,*this->delta_,*this->Adelta_,AA);
890 MVT::MvTransMv(ONE,*this->delta_,*this->Bdelta_,BB);
892 this->om_->stream(
Debug) <<
"AA: " << std::endl <<
printMat(AA) << std::endl;;
893 this->om_->stream(
Debug) <<
"BB: " << std::endl <<
printMat(BB) << std::endl;;
895 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
896 TimeMonitor lcltimer( *this->timerDS_ );
898 ret = Utils::directSolver(AA.
numRows(),AA,Teuchos::rcpFromRef(BB),S,newtheta,rank,1);
900 this->om_->stream(
Debug) <<
"S: " << std::endl <<
printMat(S) << std::endl;;
901 TEUCHOS_TEST_FOR_EXCEPTION(ret != 0,std::logic_error,
"Anasazi::IRTR::iterate(): failure solving projected eigenproblem after retraction. ret == " << ret);
908 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
909 TimeMonitor lcltimer( *this->timerSort_ );
911 std::vector<int> order(newtheta.size());
913 this->sm_->sort(newtheta, Teuchos::rcpFromRef(order), -1);
915 Utils::permuteVectors(order,S);
919 std::copy(newtheta.begin(), newtheta.begin()+this->blockSize_, this->theta_.begin());
922 this->fx_ = std::accumulate(this->theta_.begin(),this->theta_.end(),ZERO);
926 if (this->om_->isVerbosity(
Debug ) ) {
937 MagnitudeType rhonum, rhoden, mxeta;
938 std::vector<MagnitudeType> eBe(this->blockSize_);
942 rhonum = oldfx - this->fx_;
947 for (
int i=0; i<this->blockSize_; ++i) {
948 rhoden += eBe[i]*oldtheta[i];
950 mxeta = oldfx - rhoden;
951 this->rho_ = rhonum / rhoden;
952 this->om_->stream(
Debug)
953 <<
" >> old f(x) is : " << oldfx << endl
954 <<
" >> new f(x) is : " << this->fx_ << endl
955 <<
" >> m_x(eta) is : " << mxeta << endl
956 <<
" >> rhonum is : " << rhonum << endl
957 <<
" >> rhoden is : " << rhoden << endl
958 <<
" >> rho is : " << this->rho_ << endl;
960 for (
int j=0; j<this->blockSize_; ++j) {
961 this->om_->stream(
Debug)
962 <<
" >> rho[" << j <<
"] : " << 1.0/(1.0+eBe[j]) << endl;
973 this->X_ = Teuchos::null;
974 this->BX_ = Teuchos::null;
975 std::vector<int> ind(this->blockSize_);
976 for (
int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
978 X = MVT::CloneViewNonConst(*this->V_,ind);
980 BX = MVT::CloneViewNonConst(*this->BV_,ind);
982 if (useSA_ ==
false) {
986 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
987 TimeMonitor lcltimer( *this->timerLocalUpdate_ );
989 MVT::MvTimesMatAddMv(ONE,*this->delta_,S,ZERO,*X);
990 MVT::MvTimesMatAddMv(ONE,*this->Adelta_,S,ZERO,*this->AX_);
992 MVT::MvTimesMatAddMv(ONE,*this->Bdelta_,S,ZERO,*BX);
1001 Se(
Teuchos::View,S,this->blockSize_,this->blockSize_,this->blockSize_,0);
1002 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1003 TimeMonitor lcltimer( *this->timerLocalUpdate_ );
1006 MVT::MvTimesMatAddMv(ONE,* X ,Sx,ZERO,*this->delta_);
1007 MVT::MvTimesMatAddMv(ONE,*this->eta_,Se,ONE ,*this->delta_);
1008 MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*X);
1010 MVT::MvTimesMatAddMv(ONE,*this->AX_ ,Sx,ZERO,*this->delta_);
1011 MVT::MvTimesMatAddMv(ONE,*this->Aeta_,Se,ONE ,*this->delta_);
1012 MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*this->AX_);
1013 if (this->hasBOp_) {
1015 MVT::MvTimesMatAddMv(ONE,* BX ,Sx,ZERO,*this->delta_);
1016 MVT::MvTimesMatAddMv(ONE,*this->Beta_,Se,ONE ,*this->delta_);
1017 MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*BX);
1023 this->X_ = MVT::CloneView(static_cast<const MV&>(*this->V_ ),ind);
1024 this->BX_ = MVT::CloneView(static_cast<const MV&>(*this->BV_),ind);
1030 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1031 TimeMonitor lcltimer( *this->timerCompRes_ );
1033 MVT::MvAddMv( ONE, *this->BX_, ZERO, *this->BX_, *this->R_ );
1034 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
1035 MVT::MvScale( *this->R_, theta_comp );
1036 MVT::MvAddMv( ONE, *this->AX_, -ONE, *this->R_, *this->R_ );
1040 this->Rnorms_current_ =
false;
1041 this->R2norms_current_ =
false;
1046 if (this->om_->isVerbosity(
Debug ) ) {
1053 this->om_->print(
Debug, this->accuracyCheck(chk,
"after local update") );
1059 this->om_->print(
OrthoDetails, this->accuracyCheck(chk,
"after local update") );
1069 template <
class ScalarType,
class MV,
class OP>
1074 os.setf(std::ios::scientific, std::ios::floatfield);
1077 os <<
"================================================================================" << endl;
1079 os <<
" IRTR Solver Status" << endl;
1081 os <<
"The solver is "<<(this->initialized_ ?
"initialized." :
"not initialized.") << endl;
1082 os <<
"The number of iterations performed is " << this->iter_ << endl;
1083 os <<
"The current block size is " << this->blockSize_ << endl;
1084 os <<
"The number of auxiliary vectors is " << this->numAuxVecs_ << endl;
1085 os <<
"The number of operations A*x is " << this->counterAOp_ << endl;
1086 os <<
"The number of operations B*x is " << this->counterBOp_ << endl;
1087 os <<
"The number of operations B*x by the orthomanager is " << this->orthman_->getOpCounter() << endl;
1088 os <<
"The number of operations Prec*x is " << this->counterPrec_ << endl;
1089 os <<
"Parameter rho_prime is " << rho_prime_ << endl;
1090 os <<
"Inner stopping condition was " << stopReasons_[innerStop_] << endl;
1091 os <<
"Number of inner iterations was " << innerIters_ << endl;
1092 os <<
"Total number of inner iterations is " << totalInnerIters_ << endl;
1093 os <<
"f(x) is " << this->fx_ << endl;
1095 os.setf(std::ios_base::right, std::ios_base::adjustfield);
1097 if (this->initialized_) {
1099 os <<
"CURRENT EIGENVALUE ESTIMATES "<<endl;
1100 os << std::setw(20) <<
"Eigenvalue"
1101 << std::setw(20) <<
"Residual(B)"
1102 << std::setw(20) <<
"Residual(2)"
1104 os <<
"--------------------------------------------------------------------------------"<<endl;
1105 for (
int i=0; i<this->blockSize_; i++) {
1106 os << std::setw(20) << this->theta_[i];
1107 if (this->Rnorms_current_) os << std::setw(20) << this->Rnorms_[i];
1108 else os << std::setw(20) <<
"not current";
1109 if (this->R2norms_current_) os << std::setw(20) << this->R2norms_[i];
1110 else os << std::setw(20) <<
"not current";
1114 os <<
"================================================================================" << endl;
1121 #endif // ANASAZI_IRTR_HPP
This class defines the interface required by an eigensolver and status test class to compute solution...
Base class for Implicit Riemannian Trust-Region solvers.
SerialBandDenseMatrixPrinter< OrdinalType, ScalarType > printMat(const SerialBandDenseMatrix< OrdinalType, ScalarType > &obj)
Declaration of basic traits for the multivector type.
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
This class is an abstract base class for Implicit Riemannian Trust-Region based eigensolvers. The class provides the interfaces shared by the IRTR solvers (e.g., getState() and initialize()) as well as the shared implementations (e.g., inner products).
Virtual base class which defines basic traits for the operator type.
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
Pure virtual base class which describes the basic interface to the iterative eigensolver.
Anasazi's templated, static class providing utilities for the solvers.
Output managers remove the need for the eigensolver to know any information about the required output...
IRTR(const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< GenOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList ¶ms)
IRTR constructor with eigenproblem, solver utilities, and parameter list of solver options...
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
void iterate()
Impemements Eigensolver. The outer IRTR iteration. See RTRBase::iterate().
int reshape(OrdinalType numRows, OrdinalType numCols)
Types and exceptions used within Anasazi solvers and interfaces.
virtual ~IRTR()
IRTR destructor
Anasazi's templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Common interface of stopping criteria for Anasazi's solvers.
RTRRitzFailure is thrown when the RTR solver is unable to continue a call to RTRBase::iterate() due t...
OrdinalType numRows() const
void currentStatus(std::ostream &os)
Impemements Eigensolver. This method requests that the solver print out its current status to screen...