Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziIRTR.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Anasazi: Block Eigensolvers Package
4 //
5 // Copyright 2004 NTESS and the Anasazi contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
10 
11 
12 #ifndef ANASAZI_IRTR_HPP
13 #define ANASAZI_IRTR_HPP
14 
15 #include "AnasaziTypes.hpp"
16 #include "AnasaziRTRBase.hpp"
17 
18 #include "AnasaziEigensolver.hpp"
21 #include "Teuchos_ScalarTraits.hpp"
22 
23 #include "Teuchos_LAPACK.hpp"
24 #include "Teuchos_BLAS.hpp"
27 #include "Teuchos_TimeMonitor.hpp"
28 
48 // TODO: add randomization
49 // TODO: add expensive debug checking on Teuchos_Debug
50 
51 namespace Anasazi {
52 
53  template <class ScalarType, class MV, class OP>
54  class IRTR : public RTRBase<ScalarType,MV,OP> {
55  public:
56 
58 
59 
73  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
77  );
78 
80  virtual ~IRTR() {};
81 
83 
85 
86 
88  void iterate();
89 
91 
93 
94 
96  void currentStatus(std::ostream &os);
97 
99 
100  private:
101  //
102  // Convenience typedefs
103  //
104  typedef SolverUtils<ScalarType,MV,OP> Utils;
105  typedef MultiVecTraits<ScalarType,MV> MVT;
108  typedef typename SCT::magnitudeType MagnitudeType;
110  enum trRetType {
111  UNINITIALIZED = 0,
112  MAXIMUM_ITERATIONS,
113  NEGATIVE_CURVATURE,
114  EXCEEDED_TR,
115  KAPPA_CONVERGENCE,
116  THETA_CONVERGENCE
117  };
118  // these correspond to above
119  std::vector<std::string> stopReasons_;
120  //
121  // Consts
122  //
123  const MagnitudeType ZERO;
124  const MagnitudeType ONE;
125  //
126  // Internal methods
127  //
129  void solveTRSubproblem();
130  //
131  // rho_prime
132  MagnitudeType rho_prime_;
133  //
134  // norm of initial gradient: this is used for scaling
135  MagnitudeType normgradf0_;
136  //
137  // tr stopping reason
138  trRetType innerStop_;
139  //
140  // number of inner iterations
141  int innerIters_, totalInnerIters_;
142  //
143  // use 2D subspace acceleration of X+Eta to generate new iterate?
144  bool useSA_;
145  };
146 
147 
148 
149 
151  // Constructor
152  template <class ScalarType, class MV, class OP>
156  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
159  Teuchos::ParameterList &params
160  ) :
161  RTRBase<ScalarType,MV,OP>(problem,sorter,printer,tester,ortho,params,"IRTR",false),
162  ZERO(MAT::zero()),
163  ONE(MAT::one()),
164  totalInnerIters_(0)
165  {
166  // set up array of stop reasons
167  stopReasons_.push_back("n/a");
168  stopReasons_.push_back("maximum iterations");
169  stopReasons_.push_back("negative curvature");
170  stopReasons_.push_back("exceeded TR");
171  stopReasons_.push_back("kappa convergence");
172  stopReasons_.push_back("theta convergence");
173 
174  rho_prime_ = params.get("Rho Prime",0.5);
175  TEUCHOS_TEST_FOR_EXCEPTION(rho_prime_ <= 0 || rho_prime_ >= 1,std::invalid_argument,
176  "Anasazi::IRTR::constructor: rho_prime must be in (0,1).");
177 
178  useSA_ = params.get<bool>("Use SA",false);
179  }
180 
181 
183  // TR subproblem solver
184  //
185  // FINISH:
186  // define pre- and post-conditions
187  //
188  // POST:
189  // delta_,Adelta_,Bdelta_,Hdelta_ undefined
190  //
191  template <class ScalarType, class MV, class OP>
193 
194  // return one of:
195  // MAXIMUM_ITERATIONS
196  // NEGATIVE_CURVATURE
197  // EXCEEDED_TR
198  // KAPPA_CONVERGENCE
199  // THETA_CONVERGENCE
200 
201  using Teuchos::RCP;
202  using Teuchos::tuple;
203  using Teuchos::null;
204 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
205  using Teuchos::TimeMonitor;
206 #endif
207  using std::endl;
208  typedef Teuchos::RCP<const MV> PCMV;
210 
211  innerStop_ = MAXIMUM_ITERATIONS;
212 
213  const int n = MVT::GetGlobalLength(*this->eta_);
214  const int p = this->blockSize_;
215  const int d = n*p - (p*p+p)/2;
216 
217  // We have the following:
218  //
219  // X'*B*X = I
220  // X'*A*X = theta_
221  //
222  // We desire to remain in the trust-region:
223  // { eta : rho_Y(eta) \geq rho_prime }
224  // where
225  // rho_Y(eta) = 1/(1+eta'*B*eta)
226  // Therefore, the trust-region is
227  // { eta : eta'*B*eta \leq 1/rho_prime - 1 }
228  //
229  const double D2 = ONE/rho_prime_ - ONE;
230 
231  std::vector<MagnitudeType> d_Hd(p), alpha(p), beta(p), z_r(p), zold_rold(p);
232  std::vector<MagnitudeType> eBe(p), eBd(p), dBd(p), new_eBe(p);
233  MagnitudeType r0_norm;
234 
235  MVT::MvInit(*this->eta_ ,0.0);
236  MVT::MvInit(*this->Aeta_,0.0);
237  if (this->hasBOp_) {
238  MVT::MvInit(*this->Beta_,0.0);
239  }
240 
241  //
242  // R_ contains direct residuals:
243  // R_ = A X_ - B X_ diag(theta_)
244  //
245  // r0 = grad f(X) = 2 P_BX A X = 2 P_BX (A X - B X diag(theta_) = 2 proj(R_)
246  // We will do this in place.
247  // For seeking the rightmost, we want to maximize f
248  // This is the same as minimizing -f
249  // Substitute all f with -f here. In particular,
250  // grad -f(X) = -grad f(X)
251  // Hess -f(X) = -Hess f(X)
252  //
253  {
254 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
255  TimeMonitor lcltimer( *this->timerOrtho_ );
256 #endif
257  this->orthman_->projectGen(
258  *this->R_, // operating on R
259  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
260  tuple<PSDM>(null), // don't care about coeffs
261  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
262  if (this->leftMost_) {
263  MVT::MvScale(*this->R_,2.0);
264  }
265  else {
266  MVT::MvScale(*this->R_,-2.0);
267  }
268  }
269 
270  r0_norm = MAT::squareroot( RTRBase<ScalarType,MV,OP>::ginner(*this->R_) );
271  //
272  // kappa (linear) convergence
273  // theta (superlinear) convergence
274  //
275  MagnitudeType kconv = r0_norm * this->conv_kappa_;
276  // FINISH: consider inserting some scaling here
277  // MagnitudeType tconv = r0_norm * MAT::pow(r0_norm/normgradf0_,this->conv_theta_);
278  MagnitudeType tconv = MAT::pow(r0_norm,this->conv_theta_+ONE);
279  if (this->om_->isVerbosity(Debug)) {
280  this->om_->stream(Debug)
281  << " >> |r0| : " << r0_norm << endl
282  << " >> kappa conv : " << kconv << endl
283  << " >> theta conv : " << tconv << endl;
284  }
285 
286  //
287  // For Olsen preconditioning, the preconditioner is
288  // Z = P_{Prec^-1 BX, BX} Prec^-1 R
289  // for efficiency, we compute Prec^-1 BX once here for use later
290  // Otherwise, we don't need PBX
291  if (this->hasPrec_ && this->olsenPrec_)
292  {
293  std::vector<int> ind(this->blockSize_);
294  for (int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
295  Teuchos::RCP<MV> PBX = MVT::CloneViewNonConst(*this->PBV_,ind);
296  {
297 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
298  TimeMonitor prectimer( *this->timerPrec_ );
299 #endif
300  OPT::Apply(*this->Prec_,*this->BX_,*PBX);
301  this->counterPrec_ += this->blockSize_;
302  }
303  PBX = Teuchos::null;
304  }
305 
306  // Z = P_{Prec^-1 BV, BV} Prec^-1 R
307  // Prec^-1 BV in PBV
308  // or
309  // Z = P_{BV,BV} Prec^-1 R
310  if (this->hasPrec_)
311  {
312 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
313  TimeMonitor prectimer( *this->timerPrec_ );
314 #endif
315  OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
316  this->counterPrec_ += this->blockSize_;
317  // the orthogonalization time counts under Ortho and under Preconditioning
318  if (this->olsenPrec_) {
319 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
320  TimeMonitor orthtimer( *this->timerOrtho_ );
321 #endif
322  this->orthman_->projectGen(
323  *this->Z_, // operating on Z
324  tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),false, // P_{PBV,V}, B inner product, and <PBV,V>_B != I
325  tuple<PSDM>(null), // don't care about coeffs
326  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*PBV or B*Z, but do have B*V
327  }
328  else {
329 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
330  TimeMonitor orthtimer( *this->timerOrtho_ );
331 #endif
332  this->orthman_->projectGen(
333  *this->Z_, // operating on Z
334  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
335  tuple<PSDM>(null), // don't care about coeffs
336  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
337  }
338  RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,*this->Z_,z_r);
339  }
340  else {
341  RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,z_r);
342  }
343 
344  if (this->om_->isVerbosity( Debug )) {
345  // Check that gradient is B-orthogonal to X
346  typename RTRBase<ScalarType,MV,OP>::CheckList chk;
347  chk.checkBR = true;
348  if (this->hasPrec_) chk.checkZ = true;
349  this->om_->print( Debug, this->accuracyCheck(chk, "after computing gradient") );
350  }
351  else if (this->om_->isVerbosity( OrthoDetails )) {
352  // Check that gradient is B-orthogonal to X
353  typename RTRBase<ScalarType,MV,OP>::CheckList chk;
354  chk.checkBR = true;
355  if (this->hasPrec_) chk.checkZ = true;
356  this->om_->print( OrthoDetails, this->accuracyCheck(chk, "after computing gradient") );
357  }
358 
359  // delta = -z
360  MVT::MvAddMv(-ONE,*this->Z_,ZERO,*this->Z_,*this->delta_);
361 
362  if (this->om_->isVerbosity(Debug)) {
363  RCP<MV> Heta = MVT::Clone(*this->eta_,this->blockSize_);
364  // put 2*A*e - 2*B*e*theta --> He
365  MVT::MvAddMv(ONE,*this->Beta_,ZERO,*this->Beta_,*Heta);
366  std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
367  MVT::MvScale(*Heta,theta_comp);
368  if (this->leftMost_) {
369  MVT::MvAddMv( 2.0,*this->Aeta_,-2.0,*Heta,*Heta); // Heta = 2*Aeta + (-2)*Beta*theta
370  }
371  else {
372  MVT::MvAddMv(-2.0,*this->Aeta_, 2.0,*Heta,*Heta); // Heta = (-2)*Aeta + 2*Beta*theta
373  }
374  // apply projector
375  {
376 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
377  TimeMonitor lcltimer( *this->timerOrtho_ );
378 #endif
379  this->orthman_->projectGen(
380  *Heta, // operating on Heta
381  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
382  tuple<PSDM>(null), // don't care about coeffs
383  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
384  }
385 
386  std::vector<MagnitudeType> eg(this->blockSize_),
387  eHe(this->blockSize_);
388  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*this->AX_,eg);
389  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*Heta,eHe);
390  if (this->leftMost_) {
391  for (int j=0; j<this->blockSize_; ++j) {
392  eg[j] = this->theta_[j] + 2*eg[j] + .5*eHe[j];
393  }
394  }
395  else {
396  for (int j=0; j<this->blockSize_; ++j) {
397  eg[j] = -this->theta_[j] - 2*eg[j] + .5*eHe[j];
398  }
399  }
400  this->om_->stream(Debug)
401  << " Debugging checks: IRTR inner iteration " << innerIters_ << endl
402  << " >> m_X(eta) : " << std::accumulate(eg.begin(),eg.end(),0.0) << endl;
403  for (int j=0; j<this->blockSize_; ++j) {
404  this->om_->stream(Debug)
405  << " >> m_X(eta_" << j << ") : " << eg[j] << endl;
406  }
407  }
408 
411  // the inner/tCG loop
412  for (innerIters_=1; innerIters_<=d; ++innerIters_) {
413 
414  //
415  // [Hdelta,Adelta,Bdelta] = Hess*delta = 2 Proj(A*delta - B*delta*X'*A*X)
416  // X'*A*X = diag(theta), so that
417  // (B*delta)*diag(theta) can be done on the cheap
418  //
419  {
420 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
421  TimeMonitor lcltimer( *this->timerAOp_ );
422 #endif
423  OPT::Apply(*this->AOp_,*this->delta_,*this->Adelta_);
424  this->counterAOp_ += this->blockSize_;
425  }
426  if (this->hasBOp_) {
427 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
428  TimeMonitor lcltimer( *this->timerBOp_ );
429 #endif
430  OPT::Apply(*this->BOp_,*this->delta_,*this->Bdelta_);
431  this->counterBOp_ += this->blockSize_;
432  }
433  else {
434  MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*this->Bdelta_);
435  }
436  // compute <eta,B*delta> and <delta,B*delta>
437  // these will be needed below
438  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_ ,*this->Bdelta_,eBd);
439  RTRBase<ScalarType,MV,OP>::ginnersep(*this->delta_,*this->Bdelta_,dBd);
440  // put 2*A*d - 2*B*d*theta --> Hd
441  MVT::MvAddMv(ONE,*this->Bdelta_,ZERO,*this->Bdelta_,*this->Hdelta_);
442  {
443  std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
444  MVT::MvScale(*this->Hdelta_,theta_comp);
445  }
446  if (this->leftMost_) {
447  MVT::MvAddMv( 2.0,*this->Adelta_,-2.0,*this->Hdelta_,*this->Hdelta_);
448  }
449  else {
450  MVT::MvAddMv(-2.0,*this->Adelta_, 2.0,*this->Hdelta_,*this->Hdelta_);
451  }
452  // apply projector
453  {
454 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
455  TimeMonitor lcltimer( *this->timerOrtho_ );
456 #endif
457  this->orthman_->projectGen(
458  *this->Hdelta_, // operating on Hdelta
459  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
460  tuple<PSDM>(null), // don't care about coeffs
461  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
462  }
463  RTRBase<ScalarType,MV,OP>::ginnersep(*this->delta_,*this->Hdelta_,d_Hd);
464 
465 
466  // compute update step
467  for (unsigned int j=0; j<alpha.size(); ++j)
468  {
469  alpha[j] = z_r[j]/d_Hd[j];
470  }
471 
472  // compute new B-norms
473  for (unsigned int j=0; j<alpha.size(); ++j)
474  {
475  new_eBe[j] = eBe[j] + 2*alpha[j]*eBd[j] + alpha[j]*alpha[j]*dBd[j];
476  }
477 
478  if (this->om_->isVerbosity(Debug)) {
479  for (unsigned int j=0; j<alpha.size(); j++) {
480  this->om_->stream(Debug)
481  << " >> z_r[" << j << "] : " << z_r[j]
482  << " d_Hd[" << j << "] : " << d_Hd[j] << endl
483  << " >> eBe[" << j << "] : " << eBe[j]
484  << " neweBe[" << j << "] : " << new_eBe[j] << endl
485  << " >> eBd[" << j << "] : " << eBd[j]
486  << " dBd[" << j << "] : " << dBd[j] << endl;
487  }
488  }
489 
490  // check truncation criteria: negative curvature or exceeded trust-region
491  std::vector<int> trncstep;
492  trncstep.reserve(p);
493  // trncstep will contain truncated step, due to
494  // negative curvature or exceeding implicit trust-region
495  bool atleastonenegcur = false;
496  for (unsigned int j=0; j<d_Hd.size(); ++j) {
497  if (d_Hd[j] <= 0) {
498  trncstep.push_back(j);
499  atleastonenegcur = true;
500  }
501  else if (new_eBe[j] > D2) {
502  trncstep.push_back(j);
503  }
504  }
505 
506  if (!trncstep.empty())
507  {
508  // compute step to edge of trust-region, for trncstep vectors
509  if (this->om_->isVerbosity(Debug)) {
510  for (unsigned int j=0; j<trncstep.size(); ++j) {
511  this->om_->stream(Debug)
512  << " >> alpha[" << trncstep[j] << "] : " << alpha[trncstep[j]] << endl;
513  }
514  }
515  for (unsigned int j=0; j<trncstep.size(); ++j) {
516  int jj = trncstep[j];
517  alpha[jj] = ( -eBd[jj] + MAT::squareroot(eBd[jj]*eBd[jj] + dBd[jj]*(D2-eBe[jj]) ) ) / dBd[jj];
518  }
519  if (this->om_->isVerbosity(Debug)) {
520  for (unsigned int j=0; j<trncstep.size(); ++j) {
521  this->om_->stream(Debug)
522  << " >> tau[" << trncstep[j] << "] : " << alpha[trncstep[j]] << endl;
523  }
524  }
525  if (atleastonenegcur) {
526  innerStop_ = NEGATIVE_CURVATURE;
527  }
528  else {
529  innerStop_ = EXCEEDED_TR;
530  }
531  }
532 
533  // compute new eta = eta + alpha*delta
534  // we need delta*diag(alpha)
535  // do this in situ in delta_ and friends (we will note this for below)
536  // then set eta_ = eta_ + delta_
537  {
538  std::vector<ScalarType> alpha_comp(alpha.begin(),alpha.end());
539  MVT::MvScale(*this->delta_,alpha_comp);
540  MVT::MvScale(*this->Adelta_,alpha_comp);
541  if (this->hasBOp_) {
542  MVT::MvScale(*this->Bdelta_,alpha_comp);
543  }
544  MVT::MvScale(*this->Hdelta_,alpha_comp);
545  }
546  MVT::MvAddMv(ONE,*this->delta_ ,ONE,*this->eta_ ,*this->eta_);
547  MVT::MvAddMv(ONE,*this->Adelta_,ONE,*this->Aeta_,*this->Aeta_);
548  if (this->hasBOp_) {
549  MVT::MvAddMv(ONE,*this->Bdelta_,ONE,*this->Beta_,*this->Beta_);
550  }
551 
552  // store new eBe
553  eBe = new_eBe;
554 
555  //
556  // print some debugging info
557  if (this->om_->isVerbosity(Debug)) {
558  RCP<MV> Heta = MVT::Clone(*this->eta_,this->blockSize_);
559  // put 2*A*e - 2*B*e*theta --> He
560  MVT::MvAddMv(ONE,*this->Beta_,ZERO,*this->Beta_,*Heta);
561  {
562  std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
563  MVT::MvScale(*Heta,theta_comp);
564  }
565  if (this->leftMost_) {
566  MVT::MvAddMv( 2.0,*this->Aeta_,-2.0,*Heta,*Heta);
567  }
568  else {
569  MVT::MvAddMv(-2.0,*this->Aeta_, 2.0,*Heta,*Heta);
570  }
571  // apply projector
572  {
573 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
574  TimeMonitor lcltimer( *this->timerOrtho_ );
575 #endif
576  this->orthman_->projectGen(
577  *Heta, // operating on Heta
578  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
579  tuple<PSDM>(null), // don't care about coeffs
580  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
581  }
582 
583  std::vector<MagnitudeType> eg(this->blockSize_),
584  eHe(this->blockSize_);
585  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*this->AX_,eg);
586  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*Heta,eHe);
587  if (this->leftMost_) {
588  for (int j=0; j<this->blockSize_; ++j) {
589  eg[j] = this->theta_[j] + 2*eg[j] + .5*eHe[j];
590  }
591  }
592  else {
593  for (int j=0; j<this->blockSize_; ++j) {
594  eg[j] = -this->theta_[j] - 2*eg[j] + .5*eHe[j];
595  }
596  }
597  this->om_->stream(Debug)
598  << " Debugging checks: IRTR inner iteration " << innerIters_ << endl
599  << " >> m_X(eta) : " << std::accumulate(eg.begin(),eg.end(),0.0) << endl;
600  for (int j=0; j<this->blockSize_; ++j) {
601  this->om_->stream(Debug)
602  << " >> m_X(eta_" << j << ") : " << eg[j] << endl;
603  }
604  }
605 
606  //
607  // if we found negative curvature or exceeded trust-region, then quit
608  if (!trncstep.empty()) {
609  break;
610  }
611 
612  // update gradient of m
613  // R = R + Hdelta*diag(alpha)
614  // however, Hdelta_ already stores Hdelta*diag(alpha)
615  // so just add them
616  MVT::MvAddMv(ONE,*this->Hdelta_,ONE,*this->R_,*this->R_);
617  {
618  // re-tangentialize r
619 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
620  TimeMonitor lcltimer( *this->timerOrtho_ );
621 #endif
622  this->orthman_->projectGen(
623  *this->R_, // operating on R
624  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
625  tuple<PSDM>(null), // don't care about coeffs
626  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
627  }
628 
629  //
630  // check convergence
631  MagnitudeType r_norm = MAT::squareroot(RTRBase<ScalarType,MV,OP>::ginner(*this->R_,*this->R_));
632 
633  //
634  // check local convergece
635  //
636  // kappa (linear) convergence
637  // theta (superlinear) convergence
638  //
639  if (this->om_->isVerbosity(Debug)) {
640  this->om_->stream(Debug)
641  << " >> |r" << innerIters_ << "| : " << r_norm << endl;
642  }
643  if ( r_norm <= ANASAZI_MIN(tconv,kconv) ) {
644  if (tconv <= kconv) {
645  innerStop_ = THETA_CONVERGENCE;
646  }
647  else {
648  innerStop_ = KAPPA_CONVERGENCE;
649  }
650  break;
651  }
652 
653  // Z = P_{Prec^-1 BV, BV} Prec^-1 R
654  // Prec^-1 BV in PBV
655  // or
656  // Z = P_{BV,BV} Prec^-1 R
657  zold_rold = z_r;
658  if (this->hasPrec_)
659  {
660 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
661  TimeMonitor prectimer( *this->timerPrec_ );
662 #endif
663  OPT::Apply(*this->Prec_,*this->R_,*this->Z_);
664  this->counterPrec_ += this->blockSize_;
665  // the orthogonalization time counts under Ortho and under Preconditioning
666  if (this->olsenPrec_) {
667 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
668  TimeMonitor orthtimer( *this->timerOrtho_ );
669 #endif
670  this->orthman_->projectGen(
671  *this->Z_, // operating on Z
672  tuple<PCMV>(this->PBV_),tuple<PCMV>(this->V_),false, // P_{PBV,V}, B inner product, and <PBV,V>_B != I
673  tuple<PSDM>(null), // don't care about coeffs
674  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*PBV or B*Z, but do have B*V
675  }
676  else {
677 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
678  TimeMonitor orthtimer( *this->timerOrtho_ );
679 #endif
680  this->orthman_->projectGen(
681  *this->Z_, // operating on Z
682  tuple<PCMV>(this->BV_),tuple<PCMV>(this->V_),false, // P_{BV,V}, and <BV,V>_B != I
683  tuple<PSDM>(null), // don't care about coeffs
684  null,tuple<PCMV>(null), tuple<PCMV>(this->BV_)); // don't have B*BV, but do have B*V
685  }
686  RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,*this->Z_,z_r);
687  }
688  else {
689  RTRBase<ScalarType,MV,OP>::ginnersep(*this->R_,z_r);
690  }
691 
692  // compute new search direction
693  // below, we need to perform
694  // delta = -Z + delta*diag(beta)
695  // however, delta_ currently stores delta*diag(alpha)
696  // therefore, set
697  // beta_ to beta/alpha
698  // so that
699  // delta_ = delta_*diag(beta_)
700  // will in fact result in
701  // delta_ = delta_*diag(beta_)
702  // = delta*diag(alpha)*diag(beta/alpha)
703  // = delta*diag(beta)
704  // i hope this is numerically sound...
705  for (unsigned int j=0; j<beta.size(); ++j) {
706  beta[j] = z_r[j]/(zold_rold[j]*alpha[j]);
707  }
708  {
709  std::vector<ScalarType> beta_comp(beta.begin(),beta.end());
710  MVT::MvScale(*this->delta_,beta_comp);
711  }
712  MVT::MvAddMv(-ONE,*this->Z_,ONE,*this->delta_,*this->delta_);
713 
714  }
715  // end of the inner iteration loop
718  if (innerIters_ > d) innerIters_ = d;
719 
720  this->om_->stream(Debug)
721  << " >> stop reason is " << stopReasons_[innerStop_] << endl
722  << endl;
723 
724  } // end of solveTRSubproblem
725 
726 
728  // Eigensolver iterate() method
729  template <class ScalarType, class MV, class OP>
731 
732  using Teuchos::RCP;
733  using Teuchos::null;
734  using Teuchos::tuple;
735 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
736  using Teuchos::TimeMonitor;
737 #endif
738  using std::endl;
739  //typedef Teuchos::RCP<const MV> PCMV; // unused
740  //typedef Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PSDM; // unused
741 
742  //
743  // Allocate/initialize data structures
744  //
745  if (this->initialized_ == false) {
746  this->initialize();
747  }
748 
750  if (useSA_ == true) {
751  AA.reshape(2*this->blockSize_,2*this->blockSize_);
752  BB.reshape(2*this->blockSize_,2*this->blockSize_);
753  S.reshape(2*this->blockSize_,2*this->blockSize_);
754  }
755  else {
756  AA.reshape(this->blockSize_,this->blockSize_);
757  BB.reshape(this->blockSize_,this->blockSize_);
758  S.reshape(this->blockSize_,this->blockSize_);
759  }
760 
761  // set iteration details to invalid, as they don't have any meaning right now
762  innerIters_ = -1;
763  innerStop_ = UNINITIALIZED;
764 
765  // allocate temporary space
766  while (this->tester_->checkStatus(this) != Passed) {
767 
768  // Print information on current status
769  if (this->om_->isVerbosity(Debug)) {
770  this->currentStatus( this->om_->stream(Debug) );
771  }
772  else if (this->om_->isVerbosity(IterationDetails)) {
773  this->currentStatus( this->om_->stream(IterationDetails) );
774  }
775 
776  // increment iteration counter
777  this->iter_++;
778 
779  // solve the trust-region subproblem
780  solveTRSubproblem();
781  totalInnerIters_ += innerIters_;
782 
783  // perform debugging on eta et al.
784  if (this->om_->isVerbosity( Debug ) ) {
786  // this is the residual of the model, should still be in the tangent plane
787  chk.checkBR = true;
788  chk.checkEta = true;
789  chk.checkAEta = true;
790  chk.checkBEta = true;
791  this->om_->print( Debug, this->accuracyCheck(chk, "in iterate() after solveTRSubproblem()") );
792  this->om_->stream(Debug)
793  << " >> norm(Eta) : " << MAT::squareroot(RTRBase<ScalarType,MV,OP>::ginner(*this->eta_)) << endl
794  << endl;
795  }
796 
797  if (useSA_ == false) {
798  // compute the retraction of eta: R_X(eta) = X+eta
799  // we must accept, but we will work out of delta so that we can multiply back into X below
800  {
801 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
802  TimeMonitor lcltimer( *this->timerLocalUpdate_ );
803 #endif
804  MVT::MvAddMv(ONE,*this->X_,ONE,*this->eta_,*this->delta_);
805  MVT::MvAddMv(ONE,*this->AX_,ONE,*this->Aeta_,*this->Adelta_);
806  if (this->hasBOp_) {
807  MVT::MvAddMv(ONE,*this->BX_,ONE,*this->Beta_,*this->Bdelta_);
808  }
809  }
810  // perform some debugging on X+eta
811  if (this->om_->isVerbosity( Debug ) ) {
812  // X^T B X = I
813  // X^T B Eta = 0
814  // (X+Eta)^T B (X+Eta) = I + Eta^T B Eta
815  Teuchos::SerialDenseMatrix<int,ScalarType> XE(this->blockSize_,this->blockSize_),
816  E(this->blockSize_,this->blockSize_);
817  MVT::MvTransMv(ONE,*this->delta_,*this->Bdelta_,XE);
818  MVT::MvTransMv(ONE,*this->eta_,*this->Beta_,E);
819  this->om_->stream(Debug)
820  << " >> Error in AX+AEta == A(X+Eta) : " << Utils::errorEquality(*this->delta_,*this->Adelta_,this->AOp_) << endl
821  << " >> Error in BX+BEta == B(X+Eta) : " << Utils::errorEquality(*this->delta_,*this->Bdelta_,this->BOp_) << endl
822  << " >> norm( (X+Eta)^T B (X+Eta) ) : " << XE.normFrobenius() << endl
823  << " >> norm( Eta^T B Eta ) : " << E.normFrobenius() << endl
824  << endl;
825  }
826  }
827 
828  //
829  // perform rayleigh-ritz for X+eta or [X,eta] according to useSA_
830  // save an old copy of f(X) for rho analysis below
831  //
832  MagnitudeType oldfx = this->fx_;
833  std::vector<MagnitudeType> oldtheta(this->theta_), newtheta(AA.numRows());
834  int rank, ret;
835  rank = AA.numRows();
836  if (useSA_ == true) {
837 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
838  TimeMonitor lcltimer( *this->timerLocalProj_ );
839 #endif
840  Teuchos::SerialDenseMatrix<int,ScalarType> AA11(Teuchos::View,AA,this->blockSize_,this->blockSize_,0,0),
841  AA12(Teuchos::View,AA,this->blockSize_,this->blockSize_,0,this->blockSize_),
842  AA22(Teuchos::View,AA,this->blockSize_,this->blockSize_,this->blockSize_,this->blockSize_);
843  Teuchos::SerialDenseMatrix<int,ScalarType> BB11(Teuchos::View,BB,this->blockSize_,this->blockSize_,0,0),
844  BB12(Teuchos::View,BB,this->blockSize_,this->blockSize_,0,this->blockSize_),
845  BB22(Teuchos::View,BB,this->blockSize_,this->blockSize_,this->blockSize_,this->blockSize_);
846  MVT::MvTransMv(ONE,*this->X_ ,*this->AX_ ,AA11);
847  MVT::MvTransMv(ONE,*this->X_ ,*this->Aeta_,AA12);
848  MVT::MvTransMv(ONE,*this->eta_,*this->Aeta_,AA22);
849  MVT::MvTransMv(ONE,*this->X_ ,*this->BX_ ,BB11);
850  MVT::MvTransMv(ONE,*this->X_ ,*this->Beta_,BB12);
851  MVT::MvTransMv(ONE,*this->eta_,*this->Beta_,BB22);
852  }
853  else {
854 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
855  TimeMonitor lcltimer( *this->timerLocalProj_ );
856 #endif
857  MVT::MvTransMv(ONE,*this->delta_,*this->Adelta_,AA);
858  MVT::MvTransMv(ONE,*this->delta_,*this->Bdelta_,BB);
859  }
860  this->om_->stream(Debug) << "AA: " << std::endl << printMat(AA) << std::endl;;
861  this->om_->stream(Debug) << "BB: " << std::endl << printMat(BB) << std::endl;;
862  {
863 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
864  TimeMonitor lcltimer( *this->timerDS_ );
865 #endif
866  ret = Utils::directSolver(AA.numRows(),AA,Teuchos::rcpFromRef(BB),S,newtheta,rank,1);
867  }
868  this->om_->stream(Debug) << "S: " << std::endl << printMat(S) << std::endl;;
869  TEUCHOS_TEST_FOR_EXCEPTION(ret != 0,std::logic_error,"Anasazi::IRTR::iterate(): failure solving projected eigenproblem after retraction. ret == " << ret);
870  TEUCHOS_TEST_FOR_EXCEPTION(rank != AA.numRows(),RTRRitzFailure,"Anasazi::IRTR::iterate(): retracted iterate failed in Ritz analysis. rank == " << rank);
871 
872  //
873  // order the projected ritz values and vectors
874  // this ensures that the ritz vectors produced below are ordered
875  {
876 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
877  TimeMonitor lcltimer( *this->timerSort_ );
878 #endif
879  std::vector<int> order(newtheta.size());
880  // sort the values in newtheta
881  this->sm_->sort(newtheta, Teuchos::rcpFromRef(order), -1); // don't catch exception
882  // apply the same ordering to the primitive ritz vectors
883  Utils::permuteVectors(order,S);
884  }
885  //
886  // save the first blockSize values into this->theta_
887  std::copy(newtheta.begin(), newtheta.begin()+this->blockSize_, this->theta_.begin());
888  //
889  // update f(x)
890  this->fx_ = std::accumulate(this->theta_.begin(),this->theta_.end(),ZERO);
891 
892  //
893  // if debugging, do rho analysis before overwriting X,AX,BX
894  if (this->om_->isVerbosity( Debug ) ) {
895  //
896  // compute rho
897  // f(X) - f(newX) f(X) - f(newX)
898  // rho = ---------------- = ---------------------------
899  // m(0) - m(eta) -<2AX,eta> - .5*<Heta,eta>
900  //
901  // f(X) - f(newX)
902  // = ---------------------------------------
903  // -<2AX,eta> - <eta,Aeta> + <eta,Beta XAX>
904  //
905  MagnitudeType rhonum, rhoden, mxeta;
906  std::vector<MagnitudeType> eBe(this->blockSize_);
907  RTRBase<ScalarType,MV,OP>::ginnersep(*this->eta_,*this->Beta_,eBe);
908  //
909  // compute rhonum
910  rhonum = oldfx - this->fx_;
911  //
912  // compute rhoden
913  rhoden = -2.0*RTRBase<ScalarType,MV,OP>::ginner(*this->AX_ ,*this->eta_)
914  -RTRBase<ScalarType,MV,OP>::ginner(*this->Aeta_,*this->eta_);
915  for (int i=0; i<this->blockSize_; ++i) {
916  rhoden += eBe[i]*oldtheta[i];
917  }
918  mxeta = oldfx - rhoden;
919  this->rho_ = rhonum / rhoden;
920  this->om_->stream(Debug)
921  << " >> old f(x) is : " << oldfx << endl
922  << " >> new f(x) is : " << this->fx_ << endl
923  << " >> m_x(eta) is : " << mxeta << endl
924  << " >> rhonum is : " << rhonum << endl
925  << " >> rhoden is : " << rhoden << endl
926  << " >> rho is : " << this->rho_ << endl;
927  // compute individual rho
928  for (int j=0; j<this->blockSize_; ++j) {
929  this->om_->stream(Debug)
930  << " >> rho[" << j << "] : " << 1.0/(1.0+eBe[j]) << endl;
931  }
932  }
933 
934  // form new X as Ritz vectors, using the primitive Ritz vectors in S and
935  // either [X eta] or X+eta
936  // we will clear the const views of X,BX into V,BV and
937  // work from non-const temporary views
938  {
939  // release const views to X, BX
940  // get non-const views
941  this->X_ = Teuchos::null;
942  this->BX_ = Teuchos::null;
943  std::vector<int> ind(this->blockSize_);
944  for (int i=0; i<this->blockSize_; ++i) ind[i] = this->numAuxVecs_+i;
945  Teuchos::RCP<MV> X, BX;
946  X = MVT::CloneViewNonConst(*this->V_,ind);
947  if (this->hasBOp_) {
948  BX = MVT::CloneViewNonConst(*this->BV_,ind);
949  }
950  if (useSA_ == false) {
951  // multiply delta=(X+eta),Adelta=...,Bdelta=...
952  // by primitive Ritz vectors back into X,AX,BX
953  // compute ritz vectors, A,B products into X,AX,BX
954 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
955  TimeMonitor lcltimer( *this->timerLocalUpdate_ );
956 #endif
957  MVT::MvTimesMatAddMv(ONE,*this->delta_,S,ZERO,*X);
958  MVT::MvTimesMatAddMv(ONE,*this->Adelta_,S,ZERO,*this->AX_);
959  if (this->hasBOp_) {
960  MVT::MvTimesMatAddMv(ONE,*this->Bdelta_,S,ZERO,*BX);
961  }
962  }
963  else {
964  // compute ritz vectors, A,B products into X,AX,BX
965  // currently, X in X and eta in eta
966  // compute each result into delta, then copy to appropriate place
967  // decompose S into [Sx;Se]
968  Teuchos::SerialDenseMatrix<int,ScalarType> Sx(Teuchos::View,S,this->blockSize_,this->blockSize_,0,0),
969  Se(Teuchos::View,S,this->blockSize_,this->blockSize_,this->blockSize_,0);
970 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
971  TimeMonitor lcltimer( *this->timerLocalUpdate_ );
972 #endif
973  // X = [X eta] S = X*Sx + eta*Se
974  MVT::MvTimesMatAddMv(ONE,* X ,Sx,ZERO,*this->delta_);
975  MVT::MvTimesMatAddMv(ONE,*this->eta_,Se,ONE ,*this->delta_);
976  MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*X);
977  // AX = [AX Aeta] S = AX*Sx + Aeta*Se
978  MVT::MvTimesMatAddMv(ONE,*this->AX_ ,Sx,ZERO,*this->delta_);
979  MVT::MvTimesMatAddMv(ONE,*this->Aeta_,Se,ONE ,*this->delta_);
980  MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*this->AX_);
981  if (this->hasBOp_) {
982  // BX = [BX Beta] S = BX*Sx + Beta*Se
983  MVT::MvTimesMatAddMv(ONE,* BX ,Sx,ZERO,*this->delta_);
984  MVT::MvTimesMatAddMv(ONE,*this->Beta_,Se,ONE ,*this->delta_);
985  MVT::MvAddMv(ONE,*this->delta_,ZERO,*this->delta_,*BX);
986  }
987  }
988  // clear non-const views, restore const views
989  X = Teuchos::null;
990  BX = Teuchos::null;
991  this->X_ = MVT::CloneView(static_cast<const MV&>(*this->V_ ),ind);
992  this->BX_ = MVT::CloneView(static_cast<const MV&>(*this->BV_),ind);
993  }
994 
995  //
996  // update residual, R = AX - BX*theta
997  {
998 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
999  TimeMonitor lcltimer( *this->timerCompRes_ );
1000 #endif
1001  MVT::MvAddMv( ONE, *this->BX_, ZERO, *this->BX_, *this->R_ );
1002  std::vector<ScalarType> theta_comp(this->theta_.begin(),this->theta_.end());
1003  MVT::MvScale( *this->R_, theta_comp );
1004  MVT::MvAddMv( ONE, *this->AX_, -ONE, *this->R_, *this->R_ );
1005  }
1006  //
1007  // R has been updated; mark the norms as out-of-date
1008  this->Rnorms_current_ = false;
1009  this->R2norms_current_ = false;
1010 
1011 
1012  //
1013  // When required, monitor some orthogonalities
1014  if (this->om_->isVerbosity( Debug ) ) {
1015  // Check almost everything here
1017  chk.checkX = true;
1018  chk.checkAX = true;
1019  chk.checkBX = true;
1020  chk.checkR = true;
1021  this->om_->print( Debug, this->accuracyCheck(chk, "after local update") );
1022  }
1023  else if (this->om_->isVerbosity( OrthoDetails )) {
1025  chk.checkX = true;
1026  chk.checkR = true;
1027  this->om_->print( OrthoDetails, this->accuracyCheck(chk, "after local update") );
1028  }
1029 
1030  } // end while (statusTest == false)
1031 
1032  } // end of iterate()
1033 
1034 
1036  // Print the current status of the solver
1037  template <class ScalarType, class MV, class OP>
1038  void
1040  {
1041  using std::endl;
1042  os.setf(std::ios::scientific, std::ios::floatfield);
1043  os.precision(6);
1044  os <<endl;
1045  os <<"================================================================================" << endl;
1046  os << endl;
1047  os <<" IRTR Solver Status" << endl;
1048  os << endl;
1049  os <<"The solver is "<<(this->initialized_ ? "initialized." : "not initialized.") << endl;
1050  os <<"The number of iterations performed is " << this->iter_ << endl;
1051  os <<"The current block size is " << this->blockSize_ << endl;
1052  os <<"The number of auxiliary vectors is " << this->numAuxVecs_ << endl;
1053  os <<"The number of operations A*x is " << this->counterAOp_ << endl;
1054  os <<"The number of operations B*x is " << this->counterBOp_ << endl;
1055  os <<"The number of operations B*x by the orthomanager is " << this->orthman_->getOpCounter() << endl;
1056  os <<"The number of operations Prec*x is " << this->counterPrec_ << endl;
1057  os <<"Parameter rho_prime is " << rho_prime_ << endl;
1058  os <<"Inner stopping condition was " << stopReasons_[innerStop_] << endl;
1059  os <<"Number of inner iterations was " << innerIters_ << endl;
1060  os <<"Total number of inner iterations is " << totalInnerIters_ << endl;
1061  os <<"f(x) is " << this->fx_ << endl;
1062 
1063  os.setf(std::ios_base::right, std::ios_base::adjustfield);
1064 
1065  if (this->initialized_) {
1066  os << endl;
1067  os <<"CURRENT EIGENVALUE ESTIMATES "<<endl;
1068  os << std::setw(20) << "Eigenvalue"
1069  << std::setw(20) << "Residual(B)"
1070  << std::setw(20) << "Residual(2)"
1071  << endl;
1072  os <<"--------------------------------------------------------------------------------"<<endl;
1073  for (int i=0; i<this->blockSize_; i++) {
1074  os << std::setw(20) << this->theta_[i];
1075  if (this->Rnorms_current_) os << std::setw(20) << this->Rnorms_[i];
1076  else os << std::setw(20) << "not current";
1077  if (this->R2norms_current_) os << std::setw(20) << this->R2norms_[i];
1078  else os << std::setw(20) << "not current";
1079  os << endl;
1080  }
1081  }
1082  os <<"================================================================================" << endl;
1083  os << endl;
1084  }
1085 
1086 
1087 } // end Anasazi namespace
1088 
1089 #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&#39;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 &params)
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
Definition: AnasaziIRTR.hpp:80
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Common interface of stopping criteria for Anasazi&#39;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...