14#ifndef ANASAZI_MINRES_HPP
15#define ANASAZI_MINRES_HPP
18#include "Teuchos_TimeMonitor.hpp"
23template <
class Scalar,
class MV,
class OP>
24class PseudoBlockMinres
33 PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec = Teuchos::null);
36 void setTol(
const std::vector<Scalar>& tolerances) { tolerances_ = tolerances; };
37 void setMaxIter(
const int maxIt) { maxIt_ = maxIt; };
40 void setSol(Teuchos::RCP<MV> X) { X_ = X; };
41 void setRHS(Teuchos::RCP<const MV> B) { B_ = B; };
48 Teuchos::RCP<OP> Prec_;
50 Teuchos::RCP<const MV> B_;
51 std::vector<Scalar> tolerances_;
54 void symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r);
56#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
57 Teuchos::RCP<Teuchos::Time> AddTime_, ApplyOpTime_, ApplyPrecTime_, AssignTime_, DotTime_, LockTime_, NormTime_, ScaleTime_, TotalTime_;
63template <
class Scalar,
class MV,
class OP>
64PseudoBlockMinres<Scalar,MV,OP>::PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec) :
65 ONE(Teuchos::ScalarTraits<Scalar>::one()),
66 ZERO(Teuchos::ScalarTraits<Scalar>::zero())
68#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
69 AddTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Add Multivectors");
70 ApplyOpTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Apply Operator");
71 ApplyPrecTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Apply Preconditioner");
72 AssignTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Assignment (no locking)");
73 DotTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Compute Dot Product");
74 LockTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Lock Converged Vectors");
75 NormTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Compute Norm");
76 ScaleTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Scale Multivector");
77 TotalTime_ = Teuchos::TimeMonitor::getNewTimer(
"Anasazi: PseudoBlockMinres::Total Time");
87template <
class Scalar,
class MV,
class OP>
88void PseudoBlockMinres<Scalar,MV,OP>::solve()
90 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
91 Teuchos::TimeMonitor outertimer( *TotalTime_ );
95 int ncols = MVT::GetNumberVecs(*B_);
99 std::vector<Scalar> alpha(ncols), beta, beta1(ncols), phibar, oldBeta(ncols,ZERO), epsln(ncols,ZERO), cs(ncols,-ONE), sn(ncols,ZERO), dbar(ncols,ZERO), oldeps(ncols), delta(ncols), gbar(ncols), phi(ncols), gamma(ncols), tmpvec(ncols);
100 std::vector<int> indicesToRemove, newlyConverged, unconvergedIndices(ncols);
101 Teuchos::RCP<MV> V, Y, R1, R2, W, W1, W2, locX, scaleHelper, swapHelper;
104 V = MVT::Clone(*B_, ncols);
105 Y = MVT::Clone(*B_, ncols);
106 R1 = MVT::Clone(*B_, ncols);
107 R2 = MVT::Clone(*B_, ncols);
108 W = MVT::Clone(*B_, ncols);
109 W1 = MVT::Clone(*B_, ncols);
110 W2 = MVT::Clone(*B_, ncols);
111 scaleHelper = MVT::Clone(*B_, ncols);
114 indicesToRemove.reserve(ncols);
115 newlyConverged.reserve(ncols);
118 for(
int i=0; i<ncols; i++)
119 unconvergedIndices[i] = i;
123 locX = MVT::CloneCopy(*X_);
128 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
129 Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ );
131 OPT::Apply(*A_,*locX,*R1);
134 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
135 Teuchos::TimeMonitor lcltimer( *AddTime_ );
137 MVT::MvAddMv(ONE, *B_, -ONE, *R1, *R1);
142 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
143 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
145 MVT::Assign(*R1,*R2);
153 if(Prec_ != Teuchos::null)
155 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
156 Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ );
159 OPT::Apply(*Prec_,*R1,*Y);
163 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
164 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
171 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
172 Teuchos::TimeMonitor lcltimer( *DotTime_ );
174 MVT::MvDot(*Y,*R1, beta1);
176 for(
size_t i=0; i<beta1.size(); i++)
177 beta1[i] = sqrt(beta1[i]);
188 for(
int iter=1; iter <= maxIt_; iter++)
191 indicesToRemove.clear();
192 for(
int i=0; i<ncols; i++)
194 Scalar relres = phibar[i]/beta1[i];
199 if(relres < tolerances_[i])
201 indicesToRemove.push_back(i);
206 newNumConverged = indicesToRemove.size();
207 if(newNumConverged > 0)
209 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
210 Teuchos::TimeMonitor lcltimer( *LockTime_ );
214 newlyConverged.resize(newNumConverged);
215 for(
int i=0; i<newNumConverged; i++)
216 newlyConverged[i] = unconvergedIndices[indicesToRemove[i]];
218 Teuchos::RCP<MV> helperLocX = MVT::CloneCopy(*locX,indicesToRemove);
220 MVT::SetBlock(*helperLocX,newlyConverged,*X_);
223 if(newNumConverged == ncols)
227 std::vector<int> helperVec(ncols - newNumConverged);
229 std::set_difference(unconvergedIndices.begin(), unconvergedIndices.end(), newlyConverged.begin(), newlyConverged.end(), helperVec.begin());
230 unconvergedIndices = helperVec;
233 std::vector<int> thingsToKeep(ncols - newNumConverged);
234 helperVec.resize(ncols);
235 for(
int i=0; i<ncols; i++)
237 ncols = ncols - newNumConverged;
239 std::set_difference(helperVec.begin(), helperVec.end(), indicesToRemove.begin(), indicesToRemove.end(), thingsToKeep.begin());
242 Teuchos::RCP<MV> helperMV;
243 helperMV = MVT::CloneCopy(*V,thingsToKeep);
245 helperMV = MVT::CloneCopy(*Y,thingsToKeep);
247 helperMV = MVT::CloneCopy(*R1,thingsToKeep);
249 helperMV = MVT::CloneCopy(*R2,thingsToKeep);
251 helperMV = MVT::CloneCopy(*W,thingsToKeep);
253 helperMV = MVT::CloneCopy(*W1,thingsToKeep);
255 helperMV = MVT::CloneCopy(*W2,thingsToKeep);
257 helperMV = MVT::CloneCopy(*locX,thingsToKeep);
259 scaleHelper = MVT::Clone(*V,ncols);
263 oldeps.resize(ncols);
268 tmpvec.resize(ncols);
269 std::vector<Scalar> helperVecS(ncols);
270 for(
int i=0; i<ncols; i++)
271 helperVecS[i] = beta[thingsToKeep[i]];
273 for(
int i=0; i<ncols; i++)
274 helperVecS[i] = beta1[thingsToKeep[i]];
276 for(
int i=0; i<ncols; i++)
277 helperVecS[i] = phibar[thingsToKeep[i]];
279 for(
int i=0; i<ncols; i++)
280 helperVecS[i] = oldBeta[thingsToKeep[i]];
281 oldBeta = helperVecS;
282 for(
int i=0; i<ncols; i++)
283 helperVecS[i] = epsln[thingsToKeep[i]];
285 for(
int i=0; i<ncols; i++)
286 helperVecS[i] = cs[thingsToKeep[i]];
288 for(
int i=0; i<ncols; i++)
289 helperVecS[i] = sn[thingsToKeep[i]];
291 for(
int i=0; i<ncols; i++)
292 helperVecS[i] = dbar[thingsToKeep[i]];
296 A_->removeIndices(indicesToRemove);
302 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
303 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
307 for(
int i=0; i<ncols; i++)
308 tmpvec[i] = ONE/beta[i];
310 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
311 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
313 MVT::MvScale (*V, tmpvec);
319 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
320 Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ );
322 OPT::Apply(*A_, *V, *Y);
328 for(
int i=0; i<ncols; i++)
329 tmpvec[i] = beta[i]/oldBeta[i];
331 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
332 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
334 MVT::Assign(*R1,*scaleHelper);
337 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
338 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
340 MVT::MvScale(*scaleHelper,tmpvec);
343 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
344 Teuchos::TimeMonitor lcltimer( *AddTime_ );
346 MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y);
352 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
353 Teuchos::TimeMonitor lcltimer( *DotTime_ );
355 MVT::MvDot(*V, *Y, alpha);
359 for(
int i=0; i<ncols; i++)
360 tmpvec[i] = alpha[i]/beta[i];
362 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
363 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
365 MVT::Assign(*R2,*scaleHelper);
368 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
369 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
371 MVT::MvScale(*scaleHelper,tmpvec);
374 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
375 Teuchos::TimeMonitor lcltimer( *AddTime_ );
377 MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y);
388 if(Prec_ != Teuchos::null)
390 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
391 Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ );
394 OPT::Apply(*Prec_,*R2,*Y);
398 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
399 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
408 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
409 Teuchos::TimeMonitor lcltimer( *DotTime_ );
411 MVT::MvDot(*Y,*R2, beta);
413 for(
int i=0; i<ncols; i++)
414 beta[i] = sqrt(beta[i]);
419 for(
int i=0; i<ncols; i++)
421 delta[i] = cs[i]*dbar[i] + sn[i]*alpha[i];
422 gbar[i] = sn[i]*dbar[i] - cs[i]*alpha[i];
423 epsln[i] = sn[i]*beta[i];
424 dbar[i] = - cs[i]*beta[i];
428 for(
int i=0; i<ncols; i++)
430 symOrtho(gbar[i], beta[i], &cs[i], &sn[i], &gamma[i]);
432 phi[i] = cs[i]*phibar[i];
433 phibar[i] = sn[i]*phibar[i];
445 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
446 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
448 MVT::Assign(*W1,*scaleHelper);
451 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
452 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
454 MVT::MvScale(*scaleHelper,oldeps);
457 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
458 Teuchos::TimeMonitor lcltimer( *AddTime_ );
460 MVT::MvAddMv(ONE, *V, -ONE, *scaleHelper, *W);
463 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
464 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
466 MVT::Assign(*W2,*scaleHelper);
469 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
470 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
472 MVT::MvScale(*scaleHelper,delta);
475 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
476 Teuchos::TimeMonitor lcltimer( *AddTime_ );
478 MVT::MvAddMv(ONE, *W, -ONE, *scaleHelper, *W);
480 for(
int i=0; i<ncols; i++)
481 tmpvec[i] = ONE/gamma[i];
483 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
484 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
486 MVT::MvScale(*W,tmpvec);
492 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
493 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
495 MVT::Assign(*W,*scaleHelper);
498 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
499 Teuchos::TimeMonitor lcltimer( *ScaleTime_ );
501 MVT::MvScale(*scaleHelper,phi);
504 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
505 Teuchos::TimeMonitor lcltimer( *AddTime_ );
507 MVT::MvAddMv(ONE, *locX, ONE, *scaleHelper, *locX);
513 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
514 Teuchos::TimeMonitor lcltimer( *AssignTime_ );
516 MVT::SetBlock(*locX,unconvergedIndices,*X_);
520template <
class Scalar,
class MV,
class OP>
521void PseudoBlockMinres<Scalar,MV,OP>::symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r)
523 const Scalar absA = std::abs(a);
524 const Scalar absB = std::abs(b);
525 if ( absB == ZERO ) {
532 }
else if ( absA == ZERO ) {
536 }
else if ( absB >= absA ) {
539 *s = -ONE / sqrt( ONE+tau*tau );
541 *s = ONE / sqrt( ONE+tau*tau );
547 *c = -ONE / sqrt( ONE+tau*tau );
549 *c = ONE / sqrt( ONE+tau*tau );
Anasazi header file which uses auto-configuration information to include necessary C++ headers.
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package.