00001
00002
00003
00004 #ifndef ANASAZI_SIRTR_HPP
00005 #define ANASAZI_SIRTR_HPP
00006
00007 #include "AnasaziTypes.hpp"
00008 #include "AnasaziRTRBase.hpp"
00009
00010 #include "AnasaziEigensolver.hpp"
00011 #include "AnasaziMultiVecTraits.hpp"
00012 #include "AnasaziOperatorTraits.hpp"
00013 #include "Teuchos_ScalarTraits.hpp"
00014
00015 #include "Teuchos_LAPACK.hpp"
00016 #include "Teuchos_BLAS.hpp"
00017 #include "Teuchos_SerialDenseMatrix.hpp"
00018 #include "Teuchos_ParameterList.hpp"
00019 #include "Teuchos_TimeMonitor.hpp"
00020
00040
00041
00042
00043 namespace Anasazi {
00044
00045 template <class ScalarType, class MV, class OP>
00046 class SIRTR : public RTRBase<ScalarType,MV,OP> {
00047 public:
00048
00050
00051
00063 SIRTR( const Teuchos::RCP<Eigenproblem<ScalarType,MV,OP> > &problem,
00064 const Teuchos::RCP<SortManager<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > &sorter,
00065 const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00066 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00067 const Teuchos::RCP<GenOrthoManager<ScalarType,MV,OP> > &ortho,
00068 Teuchos::ParameterList ¶ms
00069 );
00070
00072 virtual ~SIRTR() {};
00073
00075
00077
00078
00080 void iterate();
00081
00083
00085
00086
00088 void currentStatus(std::ostream &os);
00089
00091
00092 private:
00093
00094
00095
00096 typedef SolverUtils<ScalarType,MV,OP> Utils;
00097 typedef MultiVecTraits<ScalarType,MV> MVT;
00098 typedef OperatorTraits<ScalarType,MV,OP> OPT;
00099 typedef Teuchos::ScalarTraits<ScalarType> SCT;
00100 typedef typename SCT::magnitudeType MagnitudeType;
00101 typedef Teuchos::ScalarTraits<MagnitudeType> MAT;
00102 enum trRetType {
00103 UNINITIALIZED = 0,
00104 MAXIMUM_ITERATIONS,
00105 NEGATIVE_CURVATURE,
00106 EXCEEDED_TR,
00107 KAPPA_CONVERGENCE,
00108 THETA_CONVERGENCE
00109 };
00110
00111 std::vector<std::string> stopReasons_;
00112
00113
00114
00115 const MagnitudeType ZERO;
00116 const MagnitudeType ONE;
00117
00118
00119
00121 void solveTRSubproblem();
00122
00123
00124 MagnitudeType rho_prime_;
00125
00126
00127 MagnitudeType normgradf0_;
00128
00129
00130 trRetType innerStop_;
00131
00132
00133 int innerIters_, totalInnerIters_;
00134 };
00135
00136
00137
00138
00140
00141 template <class ScalarType, class MV, class OP>
00142 SIRTR<ScalarType,MV,OP>::SIRTR(
00143 const Teuchos::RCP<Eigenproblem<ScalarType,MV,OP> > &problem,
00144 const Teuchos::RCP<SortManager<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > &sorter,
00145 const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00146 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00147 const Teuchos::RCP<GenOrthoManager<ScalarType,MV,OP> > &ortho,
00148 Teuchos::ParameterList ¶ms
00149 ) :
00150 RTRBase<ScalarType,MV,OP>(problem,sorter,printer,tester,ortho,params,"SIRTR",true),
00151 ZERO(MAT::zero()),
00152 ONE(MAT::one()),
00153 totalInnerIters_(0)
00154 {
00155
00156 stopReasons_.push_back("n/a");
00157 stopReasons_.push_back("maximum iterations");
00158 stopReasons_.push_back("negative curvature");
00159 stopReasons_.push_back("exceeded TR");
00160 stopReasons_.push_back("kappa convergence");
00161 stopReasons_.push_back("theta convergence");
00162
00163 rho_prime_ = params.get("Rho Prime",0.5);
00164 TEST_FOR_EXCEPTION(rho_prime_ <= 0 || rho_prime_ >= 1,std::invalid_argument,
00165 "Anasazi::SIRTR::constructor: rho_prime must be in (0,1).");
00166 }
00167
00168
00170
00171
00172
00173
00174
00175
00176
00177
00178 template <class ScalarType, class MV, class OP>
00179 void SIRTR<ScalarType,MV,OP>::solveTRSubproblem() {
00180
00181
00182
00183
00184
00185
00186
00187
00188 using Teuchos::RCP;
00189 using Teuchos::tuple;
00190 using Teuchos::null;
00191 using Teuchos::TimeMonitor;
00192 using std::endl;
00193 typedef Teuchos::RCP<const MV> PCMV;
00194 typedef Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PSDM;
00195
00196 innerStop_ = MAXIMUM_ITERATIONS;
00197
00198 const int n = MVT::GetVecLength(*this->eta_);
00199 const int p = this->blockSize_;
00200 const int d = n*p - (p*p+p)/2;
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214 const double D2 = ONE/rho_prime_ - ONE;
00215
00216 std::vector<MagnitudeType> d_Hd(p), alpha(p), beta(p), z_r(p), zold_rold(p);
00217 std::vector<MagnitudeType> eBe(p), eBd(p), dBd(p), new_eBe(p);
00218 MagnitudeType r0_norm;
00219
00220 MVT::MvInit(*this->eta_ ,0.0);
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234 {
00235 TimeMonitor lcltimer( *this->timerOrtho_ );
00236 this->orthman_->projectGen(
00237 *this->R_,
00238 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false,
00239 tuple<PSDM>(null),
00240 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00241 if (this->leftMost_) {
00242 MVT::MvScale(*this->R_,2.0);
00243 }
00244 else {
00245 MVT::MvScale(*this->R_,-2.0);
00246 }
00247 }
00248
00249 r0_norm = MAT::squareroot( ginner(*this->R_) );
00250
00251
00252
00253
00254 MagnitudeType kconv = r0_norm * this->conv_kappa_;
00255
00256
00257 MagnitudeType tconv = MAT::pow(r0_norm,this->conv_theta_+ONE);
00258 if (this->om_->isVerbosity(Debug)) {
00259 this->om_->stream(Debug)
00260 << " >> |r0| : " << r0_norm << endl
00261 << " >> kappa conv : " << kconv << endl
00262 << " >> theta conv : " << tconv << endl;
00263 }
00264
00265
00266
00267
00268
00269
00270 if (this->hasPrec_ && this->olsenPrec_)
00271 {
00272 std::vector<int> ind(this->blockSize_);
00273 for (int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
00274 Teuchos::RCP<MV> PBX = MVT::CloneView(*this->PBV_,ind);
00275 {
00276 TimeMonitor prectimer( *this->timerPrec_ );
00277 OPT::Apply(*this->Prec_,*this->BX_,*PBX);
00278 this->counterPrec_ += this->blockSize_;
00279 }
00280 PBX = Teuchos::null;
00281 }
00282
00283
00284
00285
00286
00287 if (this->hasPrec_)
00288 {
00289 TimeMonitor prectimer( *this->timerPrec_ );
00290 OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
00291 this->counterPrec_ += this->blockSize_;
00292
00293 if (this->olsenPrec_) {
00294 TimeMonitor orthtimer( *this->timerOrtho_ );
00295 this->orthman_->projectGen(
00296 *this->Z_,
00297 tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),false,
00298 tuple<PSDM>(null),
00299 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00300 }
00301 else {
00302 TimeMonitor orthtimer( *this->timerOrtho_ );
00303 this->orthman_->projectGen(
00304 *this->Z_,
00305 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false,
00306 tuple<PSDM>(null),
00307 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00308 }
00309 ginnersep(*this->R_,*this->Z_,z_r);
00310 }
00311 else {
00312
00313 MVT::MvAddMv(ONE,*this->R_,ZERO,*this->R_,*this->Z_);
00314 ginnersep(*this->R_,z_r);
00315 }
00316
00317 if (this->om_->isVerbosity( Debug )) {
00318
00319 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
00320 chk.checkBR = true;
00321 if (this->hasPrec_) chk.checkZ = true;
00322 this->om_->print( Debug, this->accuracyCheck(chk, "after computing gradient") );
00323 }
00324 else if (this->om_->isVerbosity( OrthoDetails )) {
00325
00326 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
00327 chk.checkBR = true;
00328 if (this->hasPrec_) chk.checkZ = true;
00329 this->om_->print( OrthoDetails, this->accuracyCheck(chk, "after computing gradient") );
00330 }
00331
00332
00333 MVT::MvAddMv(-ONE,*this->Z_,ZERO,*this->Z_,*this->delta_);
00334
00335 if (this->om_->isVerbosity(Debug)) {
00336
00337
00338
00339
00340 std::vector<MagnitudeType> eAx(this->blockSize_),
00341 d_eAe(this->blockSize_),
00342 d_eBe(this->blockSize_),
00343 d_mxe(this->blockSize_);
00344
00345 {
00346 TimeMonitor lcltimer( *this->timerAOp_ );
00347 OPT::Apply(*this->AOp_,*this->X_,*this->Z_);
00348 this->counterAOp_ += this->blockSize_;
00349 }
00350 ginnersep(*this->eta_,*this->Z_,eAx);
00351
00352 {
00353 TimeMonitor lcltimer( *this->timerAOp_ );
00354 OPT::Apply(*this->AOp_,*this->eta_,*this->Z_);
00355 this->counterAOp_ += this->blockSize_;
00356 }
00357 ginnersep(*this->eta_,*this->Z_,d_eAe);
00358
00359 if (this->hasBOp_) {
00360 TimeMonitor lcltimer( *this->timerBOp_ );
00361 OPT::Apply(*this->BOp_,*this->eta_,*this->Z_);
00362 this->counterBOp_ += this->blockSize_;
00363 }
00364 else {
00365 MVT::MvAddMv(ONE,*this->eta_,ZERO,*this->eta_,*this->Z_);
00366 }
00367 ginnersep(*this->eta_,*this->Z_,d_eBe);
00368
00369
00370 if (this->leftMost_) {
00371 for (int j=0; j<this->blockSize_; ++j) {
00372 d_mxe[j] = this->theta_[j] + 2*eAx[j] + d_eAe[j] - d_eBe[j]*this->theta_[j];
00373 }
00374 }
00375 else {
00376 for (int j=0; j<this->blockSize_; ++j) {
00377 d_mxe[j] = -this->theta_[j] - 2*eAx[j] - d_eAe[j] + d_eBe[j]*this->theta_[j];
00378 }
00379 }
00380 this->om_->stream(Debug)
00381 << " Debugging checks: SIRTR inner iteration " << innerIters_ << endl
00382 << " >> m_X(eta) : " << std::accumulate(d_mxe.begin(),d_mxe.end(),0.0) << endl;
00383 for (int j=0; j<this->blockSize_; ++j) {
00384 this->om_->stream(Debug)
00385 << " >> m_X(eta_" << j << ") : " << d_mxe[j] << endl;
00386 }
00387 }
00388
00391
00392 for (innerIters_=1; innerIters_<=d; ++innerIters_) {
00393
00394
00395
00396
00397
00398
00399 {
00400 TimeMonitor lcltimer( *this->timerAOp_ );
00401 OPT::Apply(*this->AOp_,*this->delta_,*this->Z_);
00402 this->counterAOp_ += this->blockSize_;
00403 }
00404 if (this->hasBOp_) {
00405 TimeMonitor lcltimer( *this->timerBOp_ );
00406 OPT::Apply(*this->BOp_,*this->delta_,*this->Hdelta_);
00407 this->counterBOp_ += this->blockSize_;
00408 }
00409 else {
00410 MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*this->Hdelta_);
00411 }
00412
00413
00414 ginnersep(*this->eta_ ,*this->Hdelta_,eBd);
00415 ginnersep(*this->delta_,*this->Hdelta_,dBd);
00416
00417 {
00418 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
00419 MVT::MvScale(*this->Hdelta_,theta_comp);
00420 }
00421 if (this->leftMost_) {
00422 MVT::MvAddMv( 2.0,*this->Z_,-2.0,*this->Hdelta_,*this->Hdelta_);
00423 }
00424 else {
00425 MVT::MvAddMv(-2.0,*this->Z_, 2.0,*this->Hdelta_,*this->Hdelta_);
00426 }
00427
00428 {
00429 TimeMonitor lcltimer( *this->timerOrtho_ );
00430 this->orthman_->projectGen(
00431 *this->Hdelta_,
00432 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false,
00433 tuple<PSDM>(null),
00434 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00435 }
00436 ginnersep(*this->delta_,*this->Hdelta_,d_Hd);
00437
00438
00439
00440 for (unsigned int j=0; j<alpha.size(); ++j)
00441 {
00442 alpha[j] = z_r[j]/d_Hd[j];
00443 }
00444
00445
00446 for (unsigned int j=0; j<alpha.size(); ++j)
00447 {
00448 new_eBe[j] = eBe[j] + 2*alpha[j]*eBd[j] + alpha[j]*alpha[j]*dBd[j];
00449 }
00450
00451 if (this->om_->isVerbosity(Debug)) {
00452 for (unsigned int j=0; j<alpha.size(); j++) {
00453 this->om_->stream(Debug)
00454 << " >> z_r[" << j << "] : " << z_r[j]
00455 << " d_Hd[" << j << "] : " << d_Hd[j] << endl
00456 << " >> eBe[" << j << "] : " << eBe[j]
00457 << " neweBe[" << j << "] : " << new_eBe[j] << endl
00458 << " >> eBd[" << j << "] : " << eBd[j]
00459 << " dBd[" << j << "] : " << dBd[j] << endl;
00460 }
00461 }
00462
00463
00464 std::vector<int> trncstep;
00465 trncstep.reserve(p);
00466
00467
00468 bool atleastonenegcur = false;
00469 for (unsigned int j=0; j<d_Hd.size(); ++j) {
00470 if (d_Hd[j] <= 0) {
00471 trncstep.push_back(j);
00472 atleastonenegcur = true;
00473 }
00474 else if (new_eBe[j] > D2) {
00475 trncstep.push_back(j);
00476 }
00477 }
00478
00479 if (!trncstep.empty())
00480 {
00481
00482 if (this->om_->isVerbosity(Debug)) {
00483 for (unsigned int j=0; j<trncstep.size(); ++j) {
00484 this->om_->stream(Debug)
00485 << " >> alpha[" << trncstep[j] << "] : " << alpha[trncstep[j]] << endl;
00486 }
00487 }
00488 for (unsigned int j=0; j<trncstep.size(); ++j) {
00489 int jj = trncstep[j];
00490 alpha[jj] = ( -eBd[jj] + MAT::squareroot(eBd[jj]*eBd[jj] + dBd[jj]*(D2-eBe[jj]) ) ) / dBd[jj];
00491 }
00492 if (this->om_->isVerbosity(Debug)) {
00493 for (unsigned int j=0; j<trncstep.size(); ++j) {
00494 this->om_->stream(Debug)
00495 << " >> tau[" << trncstep[j] << "] : " << alpha[trncstep[j]] << endl;
00496 }
00497 }
00498 if (atleastonenegcur) {
00499 innerStop_ = NEGATIVE_CURVATURE;
00500 }
00501 else {
00502 innerStop_ = EXCEEDED_TR;
00503 }
00504 }
00505
00506
00507
00508
00509
00510 {
00511 std::vector<ScalarType> alpha_comp(alpha.begin(),alpha.end());
00512 MVT::MvScale(*this->delta_,alpha_comp);
00513 MVT::MvScale(*this->Hdelta_,alpha_comp);
00514 }
00515 MVT::MvAddMv(ONE,*this->delta_ ,ONE,*this->eta_ ,*this->eta_);
00516
00517
00518 eBe = new_eBe;
00519
00520
00521
00522 if (this->om_->isVerbosity(Debug)) {
00523
00524
00525
00526
00527 std::vector<MagnitudeType> eAx(this->blockSize_),
00528 d_eAe(this->blockSize_),
00529 d_eBe(this->blockSize_),
00530 d_mxe(this->blockSize_);
00531
00532 {
00533 TimeMonitor lcltimer( *this->timerAOp_ );
00534 OPT::Apply(*this->AOp_,*this->X_,*this->Z_);
00535 this->counterAOp_ += this->blockSize_;
00536 }
00537 ginnersep(*this->eta_,*this->Z_,eAx);
00538
00539 {
00540 TimeMonitor lcltimer( *this->timerAOp_ );
00541 OPT::Apply(*this->AOp_,*this->eta_,*this->Z_);
00542 this->counterAOp_ += this->blockSize_;
00543 }
00544 ginnersep(*this->eta_,*this->Z_,d_eAe);
00545
00546 if (this->hasBOp_) {
00547 TimeMonitor lcltimer( *this->timerBOp_ );
00548 OPT::Apply(*this->BOp_,*this->eta_,*this->Z_);
00549 this->counterBOp_ += this->blockSize_;
00550 }
00551 else {
00552 MVT::MvAddMv(ONE,*this->eta_,ZERO,*this->eta_,*this->Z_);
00553 }
00554 ginnersep(*this->eta_,*this->Z_,d_eBe);
00555
00556
00557 if (this->leftMost_) {
00558 for (int j=0; j<this->blockSize_; ++j) {
00559 d_mxe[j] = this->theta_[j] + 2*eAx[j] + d_eAe[j] - d_eBe[j]*this->theta_[j];
00560 }
00561 }
00562 else {
00563 for (int j=0; j<this->blockSize_; ++j) {
00564 d_mxe[j] = -this->theta_[j] - 2*eAx[j] - d_eAe[j] + d_eBe[j]*this->theta_[j];
00565 }
00566 }
00567 this->om_->stream(Debug)
00568 << " Debugging checks: SIRTR inner iteration " << innerIters_ << endl
00569 << " >> m_X(eta) : " << std::accumulate(d_mxe.begin(),d_mxe.end(),0.0) << endl;
00570 for (int j=0; j<this->blockSize_; ++j) {
00571 this->om_->stream(Debug)
00572 << " >> m_X(eta_" << j << ") : " << d_mxe[j] << endl;
00573 }
00574 }
00575
00576
00577
00578 if (!trncstep.empty()) {
00579 break;
00580 }
00581
00582
00583
00584
00585
00586 MVT::MvAddMv(ONE,*this->Hdelta_,ONE,*this->R_,*this->R_);
00587 {
00588
00589 TimeMonitor lcltimer( *this->timerOrtho_ );
00590 this->orthman_->projectGen(
00591 *this->R_,
00592 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false,
00593 tuple<PSDM>(null),
00594 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00595 }
00596
00597
00598
00599 MagnitudeType r_norm = MAT::squareroot(ginner(*this->R_,*this->R_));
00600
00601
00602
00603
00604
00605
00606
00607 if (this->om_->isVerbosity(Debug)) {
00608 this->om_->stream(Debug)
00609 << " >> |r" << innerIters_ << "| : " << r_norm << endl;
00610 }
00611 if ( r_norm <= ANASAZI_MIN(tconv,kconv) ) {
00612 if (tconv <= kconv) {
00613 innerStop_ = THETA_CONVERGENCE;
00614 }
00615 else {
00616 innerStop_ = KAPPA_CONVERGENCE;
00617 }
00618 break;
00619 }
00620
00621
00622
00623
00624
00625 zold_rold = z_r;
00626 if (this->hasPrec_)
00627 {
00628 TimeMonitor prectimer( *this->timerPrec_ );
00629 OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
00630 this->counterPrec_ += this->blockSize_;
00631
00632 if (this->olsenPrec_) {
00633 TimeMonitor orthtimer( *this->timerOrtho_ );
00634 this->orthman_->projectGen(
00635 *this->Z_,
00636 tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),false,
00637 tuple<PSDM>(null),
00638 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00639 }
00640 else {
00641 TimeMonitor orthtimer( *this->timerOrtho_ );
00642 this->orthman_->projectGen(
00643 *this->Z_,
00644 tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false,
00645 tuple<PSDM>(null),
00646 null,tuple<PCMV>(null), tuple<PCMV>(this->BV_));
00647 }
00648 ginnersep(*this->R_,*this->Z_,z_r);
00649 }
00650 else {
00651
00652 MVT::MvAddMv(ONE,*this->R_,ZERO,*this->R_,*this->Z_);
00653 ginnersep(*this->R_,z_r);
00654 }
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669 for (unsigned int j=0; j<beta.size(); ++j) {
00670 beta[j] = z_r[j]/(zold_rold[j]*alpha[j]);
00671 }
00672 {
00673 std::vector<ScalarType> beta_comp(beta.begin(),beta.end());
00674 MVT::MvScale(*this->delta_,beta_comp);
00675 }
00676 MVT::MvAddMv(-ONE,*this->Z_,ONE,*this->delta_,*this->delta_);
00677
00678 }
00679
00682 if (innerIters_ > d) innerIters_ = d;
00683
00684 this->om_->stream(Debug)
00685 << " >> stop reason is " << stopReasons_[innerStop_] << endl
00686 << endl;
00687
00688 }
00689
00690
00691 #define SIRTR_GET_TEMP_MV(mv,workspace) \
00692 { \
00693 TEST_FOR_EXCEPTION(workspace.size() == 0,std::logic_error,"SIRTR: Request for workspace could not be honored."); \
00694 mv = workspace.back(); \
00695 workspace.pop_back(); \
00696 }
00697
00698 #define SIRTR_RELEASE_TEMP_MV(mv,workspace) \
00699 { \
00700 workspace.push_back(mv); \
00701 mv = Teuchos::null; \
00702 }
00703
00704
00706
00707 template <class ScalarType, class MV, class OP>
00708 void SIRTR<ScalarType,MV,OP>::iterate() {
00709
00710 using Teuchos::RCP;
00711 using Teuchos::null;
00712 using Teuchos::tuple;
00713 using Teuchos::TimeMonitor;
00714 using std::endl;
00715 typedef Teuchos::RCP<const MV> PCMV;
00716 typedef Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PSDM;
00717
00718
00719
00720
00721 if (this->initialized_ == false) {
00722 this->initialize();
00723 }
00724
00725 Teuchos::SerialDenseMatrix<int,ScalarType> AA(this->blockSize_,this->blockSize_),
00726 BB(this->blockSize_,this->blockSize_),
00727 S(this->blockSize_,this->blockSize_);
00728
00729
00730
00731
00732
00733
00734 std::vector< RCP<MV> > workspace;
00735
00736
00737 workspace.reserve(7);
00738
00739
00740 innerIters_ = -1;
00741 innerStop_ = UNINITIALIZED;
00742
00743
00744 while (this->tester_->checkStatus(this) != Passed) {
00745
00746
00747 if (this->om_->isVerbosity(Debug)) {
00748 this->currentStatus( this->om_->stream(Debug) );
00749 }
00750 else if (this->om_->isVerbosity(IterationDetails)) {
00751 this->currentStatus( this->om_->stream(IterationDetails) );
00752 }
00753
00754
00755 this->iter_++;
00756
00757
00758 solveTRSubproblem();
00759 totalInnerIters_ += innerIters_;
00760
00761
00762 if (this->om_->isVerbosity( Debug ) ) {
00763 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
00764
00765 chk.checkBR = true;
00766 chk.checkEta = true;
00767 this->om_->print( Debug, this->accuracyCheck(chk, "in iterate() after solveTRSubproblem()") );
00768 }
00769
00770
00771
00772
00773
00774
00775
00776
00777
00778
00779
00780 TEST_FOR_EXCEPTION(workspace.size() != 0,std::logic_error,"SIRTR::iterate(): workspace list should be empty.");
00781 SIRTR_RELEASE_TEMP_MV(this->delta_ ,workspace);
00782 SIRTR_RELEASE_TEMP_MV(this->Hdelta_,workspace);
00783 SIRTR_RELEASE_TEMP_MV(this->R_ ,workspace);
00784 SIRTR_RELEASE_TEMP_MV(this->Z_ ,workspace);
00785
00786
00787
00788
00789 RCP<MV> XpEta;
00790 SIRTR_GET_TEMP_MV(XpEta,workspace);
00791 {
00792 TimeMonitor lcltimer( *this->timerLocalUpdate_ );
00793 MVT::MvAddMv(ONE,*this->X_,ONE,*this->eta_,*XpEta);
00794 }
00795
00796
00797
00798
00799
00800 MagnitudeType oldfx = this->fx_;
00801 int rank, ret;
00802 rank = this->blockSize_;
00803
00804
00805 RCP<MV> AXpEta;
00806 SIRTR_GET_TEMP_MV(AXpEta,workspace);
00807 {
00808 TimeMonitor lcltimer( *this->timerAOp_ );
00809 OPT::Apply(*this->AOp_,*XpEta,*AXpEta);
00810 this->counterAOp_ += this->blockSize_;
00811 }
00812
00813 {
00814 TimeMonitor lcltimer( *this->timerLocalProj_ );
00815 MVT::MvTransMv(ONE,*XpEta,*AXpEta,AA);
00816 }
00817
00818
00819 RCP<MV> BXpEta;
00820 if (this->hasBOp_) {
00821 SIRTR_GET_TEMP_MV(BXpEta,workspace);
00822 {
00823 TimeMonitor lcltimer( *this->timerBOp_ );
00824 OPT::Apply(*this->BOp_,*XpEta,*BXpEta);
00825 this->counterBOp_ += this->blockSize_;
00826 }
00827
00828 {
00829 TimeMonitor lcltimer( *this->timerLocalProj_ );
00830 MVT::MvTransMv(ONE,*XpEta,*BXpEta,BB);
00831 }
00832 }
00833 else {
00834
00835 TimeMonitor lcltimer( *this->timerLocalProj_ );
00836 MVT::MvTransMv(ONE,*XpEta,*XpEta,BB);
00837 }
00838 this->om_->stream(Debug) << "AA: " << std::endl << AA << std::endl;;
00839 this->om_->stream(Debug) << "BB: " << std::endl << BB << std::endl;;
00840
00841
00842 std::vector<MagnitudeType> oldtheta(this->theta_);
00843 {
00844 TimeMonitor lcltimer( *this->timerDS_ );
00845 ret = Utils::directSolver(this->blockSize_,AA,Teuchos::rcpFromRef(BB),S,this->theta_,rank,1);
00846 }
00847 this->om_->stream(Debug) << "S: " << std::endl << S << std::endl;;
00848 TEST_FOR_EXCEPTION(ret != 0,std::logic_error,"Anasazi::SIRTR::iterate(): failure solving projected eigenproblem after retraction. ret == " << ret << "AA: " << AA << std::endl << "BB: " << BB << std::endl);
00849 TEST_FOR_EXCEPTION(rank != this->blockSize_,RTRRitzFailure,"Anasazi::SIRTR::iterate(): retracted iterate failed in Ritz analysis. rank == " << rank);
00850
00851
00852
00853
00854 {
00855 TimeMonitor lcltimer( *this->timerSort_ );
00856 std::vector<int> order(this->blockSize_);
00857
00858 this->sm_->sort(this->theta_, Teuchos::rcpFromRef(order), this->blockSize_);
00859
00860 Utils::permuteVectors(order,S);
00861 }
00862
00863
00864 this->fx_ = std::accumulate(this->theta_.begin(),this->theta_.end(),ZERO);
00865
00866
00867
00868 RCP<MV> AX;
00869 SIRTR_GET_TEMP_MV(AX,workspace);
00870 if (this->om_->isVerbosity( Debug ) ) {
00871
00872
00873
00874
00875
00876 MagnitudeType rhonum, rhoden, mxeta;
00877
00878
00879 rhonum = oldfx - this->fx_;
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890 {
00891
00892 TimeMonitor lcltimer( *this->timerAOp_ );
00893 OPT::Apply(*this->AOp_,*this->eta_,*AX);
00894 this->counterAOp_ += this->blockSize_;
00895 }
00896
00897 rhoden = -ginner(*this->eta_,*AX);
00898 if (this->hasBOp_) {
00899
00900 TimeMonitor lcltimer( *this->timerBOp_ );
00901 OPT::Apply(*this->BOp_,*this->eta_,*AX);
00902 this->counterBOp_ += this->blockSize_;
00903 }
00904 else {
00905
00906 MVT::MvAddMv(ONE,*this->eta_,ZERO,*this->eta_,*AX);
00907 }
00908
00909 std::vector<MagnitudeType> eBe(this->blockSize_);
00910 ginnersep(*this->eta_,*AX,eBe);
00911
00912 {
00913 std::vector<ScalarType> oldtheta_complex(oldtheta.begin(),oldtheta.end());
00914 MVT::MvScale( *AX, oldtheta_complex);
00915 }
00916
00917 rhoden += ginner(*this->eta_,*AX);
00918 {
00919 TimeMonitor lcltimer( *this->timerAOp_ );
00920 OPT::Apply(*this->AOp_,*this->X_,*AX);
00921 this->counterAOp_ += this->blockSize_;
00922 }
00923
00924 rhoden += -2.0*ginner(*AX,*this->eta_);
00925
00926 mxeta = oldfx - rhoden;
00927 this->rho_ = rhonum / rhoden;
00928 this->om_->stream(Debug)
00929 << " >> old f(x) is : " << oldfx << endl
00930 << " >> new f(x) is : " << this->fx_ << endl
00931 << " >> m_x(eta) is : " << mxeta << endl
00932 << " >> rhonum is : " << rhonum << endl
00933 << " >> rhoden is : " << rhoden << endl
00934 << " >> rho is : " << this->rho_ << endl;
00935
00936 for (int j=0; j<this->blockSize_; ++j) {
00937 this->om_->stream(Debug)
00938 << " >> rho[" << j << "] : " << 1.0/(1.0+eBe[j]) << endl;
00939 }
00940 }
00941
00942
00943 {
00944
00945 this->X_ = Teuchos::null;
00946 this->BX_ = Teuchos::null;
00947
00948 std::vector<int> ind(this->blockSize_);
00949 for (int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
00950 Teuchos::RCP<MV> X, BX;
00951 X = MVT::CloneView(*this->V_,ind);
00952 if (this->hasBOp_) {
00953 BX = MVT::CloneView(*this->BV_,ind);
00954 }
00955
00956 {
00957 TimeMonitor lcltimer( *this->timerLocalUpdate_ );
00958 MVT::MvTimesMatAddMv(ONE,* XpEta,S,ZERO,*X);
00959 MVT::MvTimesMatAddMv(ONE,*AXpEta,S,ZERO,*AX);
00960 if (this->hasBOp_) {
00961 MVT::MvTimesMatAddMv(ONE,*BXpEta,S,ZERO,*BX);
00962 }
00963 }
00964
00965 X = Teuchos::null;
00966 BX = Teuchos::null;
00967 this->X_ = MVT::CloneView(static_cast<const MV&>(*this->V_ ),ind);
00968 this->BX_ = MVT::CloneView(static_cast<const MV&>(*this->BV_),ind);
00969 }
00970
00971
00972 SIRTR_RELEASE_TEMP_MV(XpEta,workspace);
00973 SIRTR_RELEASE_TEMP_MV(AXpEta,workspace);
00974 if (this->hasBOp_) {
00975 SIRTR_RELEASE_TEMP_MV(BXpEta,workspace);
00976 }
00977
00978
00979
00980
00981
00982
00983 SIRTR_GET_TEMP_MV(this->R_,workspace);
00984 {
00985 TimeMonitor lcltimer( *this->timerCompRes_ );
00986 MVT::MvAddMv( ONE, *this->BX_, ZERO, *this->BX_, *this->R_ );
00987 {
00988 std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
00989 MVT::MvScale( *this->R_, theta_comp );
00990 }
00991 MVT::MvAddMv( ONE, *AX, -ONE, *this->R_, *this->R_ );
00992 }
00993
00994
00995 this->Rnorms_current_ = false;
00996 this->R2norms_current_ = false;
00997
00998
00999
01000 SIRTR_RELEASE_TEMP_MV(AX,workspace);
01001
01002
01003 SIRTR_GET_TEMP_MV(this->delta_,workspace);
01004 SIRTR_GET_TEMP_MV(this->Hdelta_,workspace);
01005 SIRTR_GET_TEMP_MV(this->Z_,workspace);
01006
01007
01008
01009 if (this->om_->isVerbosity( Debug ) ) {
01010
01011 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
01012 chk.checkX = true;
01013 chk.checkBX = true;
01014 chk.checkR = true;
01015 this->om_->print( Debug, this->accuracyCheck(chk, "after local update") );
01016 }
01017 else if (this->om_->isVerbosity( OrthoDetails )) {
01018 typename RTRBase<ScalarType,MV,OP>::CheckList chk;
01019 chk.checkX = true;
01020 chk.checkR = true;
01021 this->om_->print( OrthoDetails, this->accuracyCheck(chk, "after local update") );
01022 }
01023
01024 }
01025
01026 }
01027
01028
01030
01031 template <class ScalarType, class MV, class OP>
01032 void
01033 SIRTR<ScalarType,MV,OP>::currentStatus(std::ostream &os)
01034 {
01035 using std::endl;
01036 os.setf(std::ios::scientific, std::ios::floatfield);
01037 os.precision(6);
01038 os <<endl;
01039 os <<"================================================================================" << endl;
01040 os << endl;
01041 os <<" SIRTR Solver Status" << endl;
01042 os << endl;
01043 os <<"The solver is "<<(this->initialized_ ? "initialized." : "not initialized.") << endl;
01044 os <<"The number of iterations performed is " << this->iter_ << endl;
01045 os <<"The current block size is " << this->blockSize_ << endl;
01046 os <<"The number of auxiliary vectors is " << this->numAuxVecs_ << endl;
01047 os <<"The number of operations A*x is " << this->counterAOp_ << endl;
01048 os <<"The number of operations B*x is " << this->counterBOp_ << endl;
01049 os <<"The number of operations B*x by the orthomanager is " << this->orthman_->getOpCounter() << endl;
01050 os <<"The number of operations Prec*x is " << this->counterPrec_ << endl;
01051 os <<"Parameter rho_prime is " << rho_prime_ << endl;
01052 os <<"Inner stopping condition was " << stopReasons_[innerStop_] << endl;
01053 os <<"Number of inner iterations was " << innerIters_ << endl;
01054 os <<"Total number of inner iterations is " << totalInnerIters_ << endl;
01055 os <<"f(x) is " << this->fx_ << endl;
01056
01057 os.setf(std::ios_base::right, std::ios_base::adjustfield);
01058
01059 if (this->initialized_) {
01060 os << endl;
01061 os <<"CURRENT EIGENVALUE ESTIMATES "<<endl;
01062 os << std::setw(20) << "Eigenvalue"
01063 << std::setw(20) << "Residual(B)"
01064 << std::setw(20) << "Residual(2)"
01065 << endl;
01066 os <<"--------------------------------------------------------------------------------"<<endl;
01067 for (int i=0; i<this->blockSize_; i++) {
01068 os << std::setw(20) << this->theta_[i];
01069 if (this->Rnorms_current_) os << std::setw(20) << this->Rnorms_[i];
01070 else os << std::setw(20) << "not current";
01071 if (this->R2norms_current_) os << std::setw(20) << this->R2norms_[i];
01072 else os << std::setw(20) << "not current";
01073 os << endl;
01074 }
01075 }
01076 os <<"================================================================================" << endl;
01077 os << endl;
01078 }
01079
01080
01081 }
01082
01083 #endif // ANASAZI_SIRTR_HPP