10#ifndef ANASAZI_TraceMinBase_SOLMGR_HPP
11#define ANASAZI_TraceMinBase_SOLMGR_HPP
29#include "AnasaziStatusTestSpecTrans.hpp"
36#include "Teuchos_TimeMonitor.hpp"
37#include "Teuchos_FancyOStream.hpp"
78template<
class ScalarType,
class MV,
class OP>
84 typedef Teuchos::ScalarTraits<ScalarType> SCT;
85 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
86 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
139 Teuchos::ParameterList &
pl );
166 return Teuchos::tuple(_timerSolve, _timerRestarting, _timerLocking);
225 int blockSize_, numBlocks_, numRestartBlocks_;
231 MagnitudeType convTol_;
236 MagnitudeType lockTol_;
237 int maxLocked_, lockQuorum_;
238 bool useLocking_, relLockTol_;
242 enum WhenToShiftType whenToShift_;
243 MagnitudeType traceThresh_, shiftTol_;
244 enum HowToShiftType howToShift_;
245 bool useMultipleShifts_, relShiftTol_, considerClusters_;
246 std::string shiftNorm_;
250 std::string ortho_, which_;
251 enum SaddleSolType saddleSolType_;
252 bool projectAllVecs_, projectLockedVecs_, computeAllRes_, useRHSR_, useHarmonic_, noSort_;
253 MagnitudeType alpha_;
266 void setParameters(Teuchos::ParameterList &
pl)
const;
268 void printParameters(std::ostream &
os)
const;
271 const RCP<
SortManager<
typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> > &
sorter,
274 Teuchos::ParameterList &
plist
286template<
class ScalarType,
class MV,
class OP>
289 Teuchos::ParameterList &
pl ) :
300 TEUCHOS_TEST_FOR_EXCEPTION(problem_->getInitVec() == Teuchos::null,std::invalid_argument,
"Problem does not contain initial vectors to clone from.");
323 int osProc =
pl.get(
"Output Processor", 0);
326 Teuchos::RCP<Teuchos::FancyOStream>
osp;
328 if (
pl.isParameter(
"Output Stream")) {
329 osp = Teuchos::getParameter<Teuchos::RCP<Teuchos::FancyOStream> >(
pl,
"Output Stream");
336 if (
pl.isParameter(
"Verbosity")) {
337 if (Teuchos::isParameterType<int>(
pl,
"Verbosity")) {
340 verbosity = (
int)Teuchos::getParameter<Anasazi::MsgType>(
pl,
"Verbosity");
349 convTol_ =
pl.get(
"Convergence Tolerance",MT::prec());
351 "Anasazi::TraceMinBaseSolMgr: \"Convergence Tolerance\" must be nonnegative.");
353 relConvTol_ =
pl.get(
"Relative Convergence Tolerance",
true);
354 strtmp =
pl.get(
"Convergence Norm",std::string(
"2"));
356 convNorm_ = RES_2NORM;
359 convNorm_ = RES_ORTH;
363 "Anasazi::TraceMinBaseSolMgr: Invalid Convergence Norm.");
368 useLocking_ =
pl.get(
"Use Locking",
true);
369 relLockTol_ =
pl.get(
"Relative Locking Tolerance",
true);
370 lockTol_ =
pl.get(
"Locking Tolerance",convTol_/10);
373 "Anasazi::TraceMinBaseSolMgr: \"Relative Convergence Tolerance\" and \"Relative Locking Tolerance\" have different values. If you set one, you should always set the other.");
376 "Anasazi::TraceMinBaseSolMgr: \"Locking Tolerance\" must be nonnegative.");
378 strtmp =
pl.get(
"Locking Norm",std::string(
"2"));
380 lockNorm_ = RES_2NORM;
383 lockNorm_ = RES_ORTH;
387 "Anasazi::TraceMinBaseSolMgr: Invalid Locking Norm.");
392 maxLocked_ =
pl.get(
"Max Locked",problem_->getNEV());
394 "Anasazi::TraceMinBaseSolMgr: \"Max Locked\" must be strictly positive.");
401 lockQuorum_ =
pl.get(
"Locking Quorum",1);
403 "Anasazi::TraceMinBaseSolMgr: \"Locking Quorum\" must be strictly positive.");
410 strtmp =
pl.get(
"When To Shift",
"Always");
413 whenToShift_ = NEVER_SHIFT;
414 else if(
strtmp ==
"After Trace Levels")
415 whenToShift_ = SHIFT_WHEN_TRACE_LEVELS;
416 else if(
strtmp ==
"Residual Becomes Small")
417 whenToShift_ = SHIFT_WHEN_RESID_SMALL;
418 else if(
strtmp ==
"Always")
419 whenToShift_ = ALWAYS_SHIFT;
422 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"When To Shift\"; valid options are \"Never\", \"After Trace Levels\", \"Residual Becomes Small\", \"Always\".");
426 traceThresh_ =
pl.get(
"Trace Threshold", 0.02);
429 "Anasazi::TraceMinBaseSolMgr: \"Trace Threshold\" must be nonnegative.");
432 shiftTol_ =
pl.get(
"Shift Tolerance", 0.1);
435 "Anasazi::TraceMinBaseSolMgr: \"Shift Tolerance\" must be nonnegative.");
438 relShiftTol_ =
pl.get(
"Relative Shift Tolerance",
true);
441 shiftNorm_ =
pl.get(
"Shift Norm",
"2");
444 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"Shift Norm\"; valid options are \"2\" and \"M\".");
446 noSort_ =
pl.get(
"No Sorting",
false);
449 strtmp =
pl.get(
"How To Choose Shift",
"Adjusted Ritz Values");
451 if(
strtmp ==
"Largest Converged")
452 howToShift_ = LARGEST_CONVERGED_SHIFT;
453 else if(
strtmp ==
"Adjusted Ritz Values")
454 howToShift_ = ADJUSTED_RITZ_SHIFT;
455 else if(
strtmp ==
"Ritz Values")
456 howToShift_ = RITZ_VALUES_SHIFT;
457 else if(
strtmp ==
"Experimental Shift")
458 howToShift_ = EXPERIMENTAL_SHIFT;
461 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"How To Choose Shift\"; valid options are \"Largest Converged\", \"Adjusted Ritz Values\", \"Ritz Values\".");
464 considerClusters_ =
pl.get(
"Consider Clusters",
true);
467 useMultipleShifts_ =
pl.get(
"Use Multiple Shifts",
true);
473 ortho_ =
pl.get(
"Orthogonalization",
"SVQB");
475 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"Orthogonalization\"; valid options are \"DGKS\", \"SVQB\", \"ICGS\".");
477 strtmp =
pl.get(
"Saddle Solver Type",
"Projected Krylov");
479 if(
strtmp ==
"Projected Krylov")
480 saddleSolType_ = PROJECTED_KRYLOV_SOLVER;
481 else if(
strtmp ==
"Schur Complement")
482 saddleSolType_ = SCHUR_COMPLEMENT_SOLVER;
483 else if(
strtmp ==
"Block Diagonal Preconditioned Minres")
484 saddleSolType_ = BD_PREC_MINRES;
485 else if(
strtmp ==
"HSS Preconditioned Gmres")
486 saddleSolType_ = HSS_PREC_GMRES;
489 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"Saddle Solver Type\"; valid options are \"Projected Krylov\", \"Schur Complement\", and \"Block Diagonal Preconditioned Minres\".");
491 projectAllVecs_ =
pl.get(
"Project All Vectors",
true);
492 projectLockedVecs_ =
pl.get(
"Project Locked Vectors",
true);
493 computeAllRes_ =
pl.get(
"Compute All Residuals",
true);
494 useRHSR_ =
pl.get(
"Use Residual as RHS",
false);
495 alpha_ =
pl.get(
"HSS: alpha", 1.0);
498 "Anasazi::TraceMinBaseSolMgr: If you want to project out the locked vectors, you should really project out ALL the vectors of X.");
501 maxKrylovIter_ =
pl.get(
"Maximum Krylov Iterations", 200);
503 "Anasazi::TraceMinBaseSolMgr: \"Maximum Krylov Iterations\" must be greater than 0.");
506 which_ =
pl.get(
"Which",
"SM");
508 "Anasazi::TraceMinBaseSolMgr: Invalid value for \"Which\"; valid options are \"SM\" and \"LM\".");
513 "Anasazi::TraceMinBaseSolMgr: It is an exceptionally bad idea to use Ritz shifts when finding the largest eigenpairs of a standard eigenvalue problem. If we don't use Ritz shifts, it may take extra iterations to converge, but we NEVER have to solve a single linear system. Using Ritz shifts forces us to solve systems of the form (I + sigma A)x=f, and it probably doesn't benefit us enough to outweigh the extra cost. We may add support for this feature in the future, but for now, please set \"When To Shift\" to \"Never\".");
515#ifdef BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
518 TEUCHOS_TEST_FOR_EXCEPTION(problem_->getPrec() != Teuchos::null && saddleSolType_ == PROJECTED_KRYLOV_SOLVER && useMultipleShifts_, std::invalid_argument,
519 "Anasazi::TraceMinBaseSolMgr: When you use the projected Krylov solver with preconditioning, the preconditioner must be projected as well. In theory, if the preconditioner is SPD, the projected preconditioner will be SPSD, but in practice, it can have small negative eigenvalues, presumably due to machine arithmetic. This means we can't use TraceMin's built-in MINRES, and we are forced to use Belos for now. When you use multiple Ritz shifts, you are essentially using a different operator to solve each linear system. Belos can't handle this right now, but we're working on a solution. For now, please set \"Use Multiple Shifts\" to false.");
526 TEUCHOS_TEST_FOR_EXCEPTION(problem_->getPrec() != Teuchos::null && saddleSolType_ == PROJECTED_KRYLOV_SOLVER, std::invalid_argument,
527 "Anasazi::TraceMinBaseSolMgr: When you use the projected Krylov solver with preconditioning, the preconditioner must be projected as well. In theory, if the preconditioner is SPD, the projected preconditioner will be SPSD, but in practice, it can have small negative eigenvalues, presumably due to machine arithmetic. This means we can't use TraceMin's built-in MINRES, and we are forced to use Belos for now. You didn't install Belos. You have three options to correct this problem:\n1. Reinstall Trilinos with Belos enabled\n2. Don't use a preconditioner\n3. Choose a different method for solving the saddle-point problem (Recommended)");
539template<
class ScalarType,
class MV,
class OP>
545 const int nev = problem_->getNEV();
549 out = Teuchos::getFancyOStream(Teuchos::rcpFromRef(printer_->stream(
Debug)));
550 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
551 *
out <<
"Entering Anasazi::TraceMinBaseSolMgr::solve()\n";
564 problem_->setOperator(problem_->getM());
566 problem_->setProblem();
574 if (globalTest_ == Teuchos::null) {
588 if (lockingTest_ == Teuchos::null) {
599 Teuchos::Array<RCP<StatusTest<ScalarType,MV,OP> > >
alltests;
602 if (debugTest_ != Teuchos::null)
alltests.push_back(debugTest_);
608 if ( printer_->isVerbosity(
Debug) ) {
618 if (ortho_==
"SVQB") {
620 }
else if (ortho_==
"DGKS") {
622 }
else if (ortho_==
"ICGS") {
630 Teuchos::ParameterList
plist;
631 setParameters(
plist);
638 RCP< const MV >
probauxvecs = problem_->getAuxVecs();
654 if (maxLocked_ > 0) {
655 lockvecs = MVT::Clone(*problem_->getInitVec(),maxLocked_);
657 std::vector<MagnitudeType>
lockvals;
702 const ScalarType ONE = SCT::one();
703 const ScalarType ZERO = SCT::zero();
708 problem_->setSolution(
sol);
714#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
715 Teuchos::TimeMonitor
slvtimer(*_timerSolve);
727 if (debugTest_ != Teuchos::null && debugTest_->getStatus() ==
Passed) {
728 throw AnasaziError(
"Anasazi::TraceMinBaseSolMgr::solve(): User-specified debug status test returned Passed.");
748#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
749 Teuchos::TimeMonitor
lcktimer(*_timerLocking);
760 "Anasazi::TraceMinBaseSolMgr::solve(): status test mistake: howMany() non-positive.");
762 "Anasazi::TraceMinBaseSolMgr::solve(): status test mistake: howMany() not consistent with whichVecs().");
765 "Anasazi::TraceMinBaseSolMgr::solve(): status test mistake: locking not deactivated.");
795 if (printer_->isVerbosity(
Debug)) {
796 printer_->print(
Debug,
"Locking vectors: ");
798 printer_->print(
Debug,
"\n");
799 printer_->print(
Debug,
"Not locking vectors: ");
801 printer_->print(
Debug,
"\n");
807 printer_->stream(
Debug) <<
"Ritz value[" <<
i <<
"] = " <<
allvals[
i].realpart << std::endl;
849 Teuchos::Array< RCP<const MV> >
aux;
856 printer_->stream(
Debug) <<
"maxLocked: " << maxLocked_ << std::endl;
861 printer_->stream(
Debug) <<
"Removed locking test\n";
864 int newdim = numRestartBlocks_*blockSize_;
919 MVT::MvRandom(*
newV);
926 for(
unsigned int i=0;
i<(
unsigned int)blockSize_;
i++)
964 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"Anasazi::TraceMinBaseSolMgr::solve(): Invalid return from tm_solver::iterate().");
969 <<
"Anasazi::TraceMinBaseSolMgr::solve() caught unexpected exception from Anasazi::TraceMinBase::iterate() at iteration " <<
tm_solver->getNumIters() << std::endl
970 <<
err.what() << std::endl
971 <<
"Anasazi::TraceMinBaseSolMgr::solve() returning Unconverged with no solutions." << std::endl;
977 if (
sol.numVecs > 0) {
978 sol.Evecs = MVT::Clone(*problem_->getInitVec(),
sol.numVecs);
980 sol.Evals.resize(
sol.numVecs);
981 std::vector<MagnitudeType>
vals(
sol.numVecs);
989 for (
unsigned int i=0;
i<
which.size();
i++) {
1037 for(
size_t i=0;
i<
vals.size();
i++)
1043 std::vector<int>
order(
sol.numVecs);
1046 for (
int i=0;
i<
sol.numVecs;
i++) {
1048 sol.Evals[
i].imagpart = MT::zero();
1051 msutils::permuteVectors(
sol.numVecs,
order,*
sol.Evecs);
1055 sol.index.resize(
sol.numVecs,0);
1065#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1067 Teuchos::TimeMonitor::summarize( printer_->stream(
TimingDetails ) );
1071 problem_->setSolution(
sol);
1072 printer_->stream(
Debug) <<
"Returning " <<
sol.numVecs <<
" eigenpairs to eigenproblem." << std::endl;
1084template <
class ScalarType,
class MV,
class OP>
1092template <
class ScalarType,
class MV,
class OP>
1099template <
class ScalarType,
class MV,
class OP>
1107template <
class ScalarType,
class MV,
class OP>
1114template <
class ScalarType,
class MV,
class OP>
1122template <
class ScalarType,
class MV,
class OP>
1126 return lockingTest_;
1129template <
class ScalarType,
class MV,
class OP>
1132 const ScalarType ONE = Teuchos::ScalarTraits<MagnitudeType>::one();
1133 const ScalarType ZERO = Teuchos::ScalarTraits<MagnitudeType>::zero();
1165 if(problem_->getOperator() != Teuchos::null)
1172 newState.KV = Teuchos::null;
1173 newState.KX = Teuchos::null;
1176 if(problem_->getM() != Teuchos::null)
1178 newState.MopV = MVT::CloneView(*oldState.MX, indToCopy);
1179 newState.MX = newState.MopV;
1183 newState.MopV = Teuchos::null;
1184 newState.MX = Teuchos::null;
1187 else if(newdim == 0)
1189 std::vector<int> blockind(blockSize_);
1190 for(
int i=0; i<blockSize_; i++)
1194 newState.V = MVT::CloneView(*oldState.X, blockind);
1195 newState.KV = MVT::CloneView(*oldState.KX, blockind);
1196 newState.R = MVT::CloneView(*oldState.R, blockind);
1197 newState.X = MVT::CloneView(*newState.V, blockind);
1198 newState.KX = MVT::CloneView(*newState.KV, blockind);
1200 if(problem_->getM() != Teuchos::null)
1202 newState.MopV = MVT::CloneView(*oldState.MX, blockind);
1203 newState.MX = MVT::CloneView(*newState.MopV, blockind);
1207 newState.MopV = Teuchos::null;
1208 newState.MX = Teuchos::null;
1214 std::vector<int> oldPart(olddim);
1215 for(
int i=0; i<olddim; i++) oldPart[i] = i;
1216 std::vector<int> newPart(newdim);
1217 for(
int i=0; i<newdim; i++) newPart[i] = olddim+i;
1220 RCP<MV> helper = MVT::Clone(*oldState.V,newState.curDim);
1221 RCP<MV> oldHelper = MVT::CloneViewNonConst(*helper,oldPart);
1222 RCP<MV> newHelper = MVT::CloneViewNonConst(*helper,newPart);
1223 RCP<const MV> viewHelper;
1226 Teuchos::SerialDenseMatrix<int,ScalarType> newRV(oldState.curDim,newdim);
1227 for(
int r=0; r<oldState.curDim; r++)
1229 for(
int c=0; c<newdim; c++)
1230 newRV(r,c) = (*oldState.RV)(r,newIndices[c]);
1234 viewHelper = MVT::CloneView(*oldState.V,fullIndices);
1235 MVT::MvTimesMatAddMv(ONE,*viewHelper,newRV,ZERO,*newHelper);
1236 viewHelper = MVT::CloneView(*oldState.X,oldIndices);
1237 MVT::Assign(*viewHelper,*oldHelper);
1238 newState.V = MVT::CloneCopy(*helper);
1241 viewHelper = MVT::CloneView(*oldState.KV,fullIndices);
1242 MVT::MvTimesMatAddMv(ONE,*viewHelper,newRV,ZERO,*newHelper);
1243 viewHelper = MVT::CloneView(*oldState.KX,oldIndices);
1244 MVT::Assign(*viewHelper,*oldHelper);
1245 newState.KV = MVT::CloneCopy(*helper);
1248 if(problem_->getM() != Teuchos::null)
1250 viewHelper = MVT::CloneView(*oldState.MopV,fullIndices);
1251 MVT::MvTimesMatAddMv(ONE,*viewHelper,newRV,ZERO,*newHelper);
1252 viewHelper = MVT::CloneView(*oldState.MX,oldIndices);
1253 MVT::Assign(*viewHelper,*oldHelper);
1254 newState.MopV = MVT::CloneCopy(*helper);
1257 newState.MopV = newState.V;
1260 std::vector<int> blockVec(blockSize_);
1261 for(
int i=0; i<blockSize_; i++) blockVec[i] = i;
1262 newState.X = MVT::CloneView(*newState.V,blockVec);
1263 newState.KX = MVT::CloneView(*newState.KV,blockVec);
1264 newState.MX = MVT::CloneView(*newState.MopV,blockVec);
1267 if(blockSize_-oldIndices.size() > 0)
1270 newPart.resize(blockSize_-oldIndices.size());
1271 helper = MVT::Clone(*oldState.V,blockSize_);
1272 oldHelper = MVT::CloneViewNonConst(*helper,oldPart);
1273 newHelper = MVT::CloneViewNonConst(*helper,newPart);
1275 RCP<MV> scaledMV = MVT::CloneCopy(*newState.MX,newPart);
1276 RCP<const MV> localKV = MVT::CloneView(*newState.KX,newPart);
1277 std::vector<ScalarType> scalarVec(blockSize_-oldIndices.size());
1278 for(
unsigned int i=0; i<(
unsigned int)blockSize_-oldIndices.size(); i++) scalarVec[i] = (*oldState.T)[newPart[i]];
1279 MVT::MvScale(*scaledMV,scalarVec);
1281 helper = MVT::Clone(*oldState.V,blockSize_);
1282 oldHelper = MVT::CloneViewNonConst(*helper,oldPart);
1283 newHelper = MVT::CloneViewNonConst(*helper,newPart);
1284 MVT::MvAddMv(ONE,*localKV,-ONE,*scaledMV,*newHelper);
1285 viewHelper = MVT::CloneView(*oldState.R,oldIndices);
1286 MVT::Assign(*viewHelper,*oldHelper);
1287 newState.R = MVT::CloneCopy(*helper);
1290 newState.R = oldState.R;
1294 newState.isOrtho =
true;
1297 RCP< std::vector<ScalarType> > helperT = rcp(
new std::vector<ScalarType>(newState.curDim) );
1298 for(
int i=0; i<newState.curDim; i++) (*helperT)[i] = (*oldState.T)[indToCopy[i]];
1299 newState.T = helperT;
1302 RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > newKK = rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(newState.curDim,newState.curDim) );
1303 for(
int i=0; i<newState.curDim; i++)
1304 (*newKK)(i,i) = (*newState.T)[i];
1305 newState.KK = newKK;
1308 RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > newRV = rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(newState.curDim,newState.curDim) );
1309 for(
int i=0; i<newState.curDim; i++)
1310 (*newRV)(i,i) = ONE;
1311 newState.RV = newRV;
1314 RCP< std::vector<ScalarType> > helperRS = rcp(
new std::vector<ScalarType>(blockSize_) );
1315 for(
int i=0; i<blockSize_; i++)
1317 if(indToCopy[i] < blockSize_)
1318 (*helperRS)[i] = (*oldState.ritzShifts)[indToCopy[i]];
1320 (*helperRS)[i] = ZERO;
1322 newState.ritzShifts = helperRS;
1326template <
class ScalarType,
class MV,
class OP>
1327void TraceMinBaseSolMgr<ScalarType,MV,OP>::setParameters(Teuchos::ParameterList &pl)
const
1329 pl.set(
"Block Size", blockSize_);
1330 pl.set(
"Num Blocks", numBlocks_);
1331 pl.set(
"Num Restart Blocks", numRestartBlocks_);
1332 pl.set(
"When To Shift", whenToShift_);
1333 pl.set(
"Trace Threshold", traceThresh_);
1334 pl.set(
"Shift Tolerance", shiftTol_);
1335 pl.set(
"Relative Shift Tolerance", relShiftTol_);
1336 pl.set(
"Shift Norm", shiftNorm_);
1337 pl.set(
"How To Choose Shift", howToShift_);
1338 pl.set(
"Consider Clusters", considerClusters_);
1339 pl.set(
"Use Multiple Shifts", useMultipleShifts_);
1340 pl.set(
"Saddle Solver Type", saddleSolType_);
1341 pl.set(
"Project All Vectors", projectAllVecs_);
1342 pl.set(
"Project Locked Vectors", projectLockedVecs_);
1343 pl.set(
"Compute All Residuals", computeAllRes_);
1344 pl.set(
"Use Residual as RHS", useRHSR_);
1345 pl.set(
"Use Harmonic Ritz Values", useHarmonic_);
1346 pl.set(
"Maximum Krylov Iterations", maxKrylovIter_);
1347 pl.set(
"HSS: alpha", alpha_);
1351template <
class ScalarType,
class MV,
class OP>
1352void TraceMinBaseSolMgr<ScalarType,MV,OP>::printParameters(std::ostream &os)
const
1355 os <<
"========================================\n";
1356 os <<
"========= TraceMin parameters ==========\n";
1357 os <<
"========================================\n";
1358 os <<
"=========== Block parameters ===========\n";
1359 os <<
"Block Size: " << blockSize_ << std::endl;
1360 os <<
"Num Blocks: " << numBlocks_ << std::endl;
1361 os <<
"Num Restart Blocks: " << numRestartBlocks_ << std::endl;
1362 os <<
"======== Convergence parameters ========\n";
1363 os <<
"Convergence Tolerance: " << convTol_ << std::endl;
1364 os <<
"Relative Convergence Tolerance: " << relConvTol_ << std::endl;
1365 os <<
"========== Locking parameters ==========\n";
1366 os <<
"Use Locking: " << useLocking_ << std::endl;
1367 os <<
"Locking Tolerance: " << lockTol_ << std::endl;
1368 os <<
"Relative Locking Tolerance: " << relLockTol_ << std::endl;
1369 os <<
"Max Locked: " << maxLocked_ << std::endl;
1370 os <<
"Locking Quorum: " << lockQuorum_ << std::endl;
1371 os <<
"========== Shifting parameters =========\n";
1372 os <<
"When To Shift: ";
1373 if(whenToShift_ == NEVER_SHIFT) os <<
"Never\n";
1374 else if(whenToShift_ == SHIFT_WHEN_TRACE_LEVELS) os <<
"After Trace Levels\n";
1375 else if(whenToShift_ == SHIFT_WHEN_RESID_SMALL) os <<
"Residual Becomes Small\n";
1376 else if(whenToShift_ == ALWAYS_SHIFT) os <<
"Always\n";
1377 os <<
"Consider Clusters: " << considerClusters_ << std::endl;
1378 os <<
"Trace Threshohld: " << traceThresh_ << std::endl;
1379 os <<
"Shift Tolerance: " << shiftTol_ << std::endl;
1380 os <<
"Relative Shift Tolerance: " << relShiftTol_ << std::endl;
1381 os <<
"How To Choose Shift: ";
1382 if(howToShift_ == LARGEST_CONVERGED_SHIFT) os <<
"Largest Converged\n";
1383 else if(howToShift_ == ADJUSTED_RITZ_SHIFT) os <<
"Adjusted Ritz Values\n";
1384 else if(howToShift_ == RITZ_VALUES_SHIFT) os <<
"Ritz Values\n";
1385 os <<
"Use Multiple Shifts: " << useMultipleShifts_ << std::endl;
1386 os <<
"=========== Other parameters ===========\n";
1387 os <<
"Orthogonalization: " << ortho_ << std::endl;
1388 os <<
"Saddle Solver Type: ";
1389 if(saddleSolType_ == PROJECTED_KRYLOV_SOLVER) os <<
"Projected Krylov\n";
1390 else if(saddleSolType_ == SCHUR_COMPLEMENT_SOLVER) os <<
"Schur Complement\n";
1391 os <<
"Project All Vectors: " << projectAllVecs_ << std::endl;
1392 os <<
"Project Locked Vectors: " << projectLockedVecs_ << std::endl;
1393 os <<
"Compute All Residuals: " << computeAllRes_ << std::endl;
1394 os <<
"========================================\n\n\n";
Basic implementation of the Anasazi::OrthoManager class.
Basic implementation of the Anasazi::SortManager class.
Anasazi header file which uses auto-configuration information to include necessary C++ headers.
Abstract base class which defines the interface required by an eigensolver and status test class to c...
Basic implementation of the Anasazi::OrthoManager class.
Abstract class definition for Anasazi Output Managers.
Abstract class definition for Anasazi output stream.
Orthogonalization manager based on the SVQB technique described in "A Block Orthogonalization Procedu...
Pure virtual base class which describes the basic interface for a solver manager.
Class which provides internal utilities for the Anasazi solvers.
Status test for forming logical combinations of other status tests.
Special StatusTest for printing status tests.
A status test for testing the norm of the eigenvectors residuals.
A status test for testing the norm of the eigenvectors residuals along with a set of auxiliary eigenv...
Abstract base class for trace minimization eigensolvers.
Types and exceptions used within Anasazi solvers and interfaces.
An exception class parent to all Anasazi exceptions.
The Anasazi::TraceMinBaseSolMgr provides an abstract base class for the TraceMin series of solver man...
const RCP< StatusTest< ScalarType, MV, OP > > & getLockingStatusTest() const
Get the status test defining locking.
Teuchos::Array< RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
const RCP< StatusTest< ScalarType, MV, OP > > & getDebugStatusTest() const
Get the status test for debugging.
ReturnType solve()
This method performs possibly repeated calls to the underlying eigensolver's iterate() routine until ...
void setDebugStatusTest(const RCP< StatusTest< ScalarType, MV, OP > > &debug)
Set the status test for debugging.
void setGlobalStatusTest(const RCP< StatusTest< ScalarType, MV, OP > > &global)
Set the status test defining global convergence.
virtual ~TraceMinBaseSolMgr()
Destructor.
TraceMinBaseSolMgr(const RCP< Eigenproblem< ScalarType, MV, OP > > &problem, Teuchos::ParameterList &pl)
Basic constructor for TraceMinBaseSolMgr.
int getNumIters() const
Get the iteration count for the most recent call to solve().
void setLockingStatusTest(const RCP< StatusTest< ScalarType, MV, OP > > &locking)
Set the status test defining locking.
const RCP< StatusTest< ScalarType, MV, OP > > & getGlobalStatusTest() const
Get the status test defining global convergence.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Return the eigenvalue problem.
Anasazi's templated virtual class for constructing an operator that can interface with the OperatorTr...
The Anasazi::SolverManager is a templated virtual base class that defines the basic interface that an...
Anasazi's templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Status test for forming logical combinations of other status tests.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package.
ReturnType
Enumerated type used to pass back information from a solver manager.
ResType
Enumerated type used to specify which residual norm used by residual norm status tests.
Output managers remove the need for the eigensolver to know any information about the required output...