Class documentation of Concepts

Loading...
Searching...
No Matches
anasaziES.hh
Go to the documentation of this file.
1
5#pragma once
6
9#include "operator.hh"
10#include <AnasaziBasicEigenproblem.hpp>
11#include <AnasaziSolverManager.hpp>
12#include <Teuchos_ParameterList.hpp>
13
14namespace concepts {
15
16template <class ScalarT>
17class AnasaziES : public eigensolver::EigenSolver<ScalarT> {
18public:
19 enum State {INIT_PHASE, CONFIGURED_PHASE, SOLVED_PHASE};
20 enum SolverType {BLOCK_KRYLOV_SCHUR, BLOCK_DAVIDSON, LOBPCG, RTR};
21
22 typedef Anasazi::MultiVec<ScalarT> MV; // abstract base
23 typedef Anasazi::Operator<ScalarT> OP; // abstract base
24 typedef AnasaziMV<ScalarT> MV_C; // concrete implementation
25 typedef AnasaziOp<ScalarT> OP_C; // concrete implementation
26 typedef Anasazi::MultiVecTraits<ScalarT, MV> MVT;
27
28 AnasaziES(SolverType solverType,
30 );
31
32 ~AnasaziES() {
33
34 }
35
36 void setM(concepts::Operator<ScalarT>& M) {
37 assert(phase_ == INIT_PHASE);
38 myProblem_->setM(Teuchos::rcp( new OP_C(Teuchos::rcp(&M, false) ) ));
39 }
40
41 void setHermitian(bool isHermite) {
42 assert(phase_ == INIT_PHASE);
43 myProblem_->setHermitian(isHermite);
44 }
45
47 void setBlockSize(int blocksize = 1) {
48 myPL_.set( "Block Size", blocksize );
49 int dim = A_->getOp().dimX();
50
51 Teuchos::RCP<MV> ivec = Teuchos::rcp( new MV_C(dim,blocksize) );
52 ivec->MvRandom();
53 myProblem_->setInitVec(ivec);
54 }
55
56 void setNEV(int nev=3) {
57 assert(phase_ == INIT_PHASE);
58 myProblem_->setNEV(nev);
59 }
60
66 bool setProblem();
67
68 int solve();
69
70 /*
71 Specify the verbosity level. Options include:
72 Anasazi::Errors
73 This option is always set
74 Anasazi::Warnings
75 Warnings (less severe than errors)
76 Anasazi::IterationDetails
77 Details at each iteration, such as the current eigenvalues
78 Anasazi::OrthoDetails
79 Details about orthogonality
80 Anasazi::TimingDetails
81 A summary of the timing info for the solve() routine
82 Anasazi::FinalSummary
83 A final summary
84 Anasazi::Debug
85 Debugging information
86 */
87 void setVerbosity(int verbosity = Anasazi::Errors) {
88 assert(phase_ == INIT_PHASE);
89 myPL_.set( "Verbosity", verbosity );
90 }
91
93 void setOrthogonalization(std::string mode = "SVQB") {
94 assert(phase_ == INIT_PHASE);
95 myPL_.set( "Orthogonalization", mode );
96 }
97
98 /*
99 Choose which eigenvalues to compute.
100 Choices are:
101 LM - target the largest magnitude [default]
102 SM - target the smallest magnitude
103 LR - target the largest real
104 SR - target the smallest real
105 LI - target the largest imaginary
106 SI - target the smallest imaginary
107 */
108 void setMode(std::string mode = "LM") {
109 assert(phase_ == INIT_PHASE);
110 myPL_.set( "Which", mode );
111 }
112
117 assert(phase_ == INIT_PHASE);
118 myPL_.set( "Num Blocks", numBlocks);
119 }
120
121 void setMaxRestarts(int maxRestarts = 100) {
122 assert(phase_ == INIT_PHASE);
123 myPL_.set( "Maximum Restarts", maxRestarts );
124 }
125
129 void setConvTol(double tolerance) {
130 assert(phase_ == INIT_PHASE);
131 myPL_.set( "Convergence Tolerance", tolerance);
132 }
133
138 void setRelConvTol(bool rct) {
139 assert(phase_ == INIT_PHASE);
140 myPL_.set( "Relative Convergence Tolerance", rct);
141 }
142
145 void setConvergenceNorm(std::string norm = "2") {
146 assert(phase_ == INIT_PHASE);
147 myPL_.set( "Convergence Norm", norm );
148 }
149
153 void setMaxLocked(int val) {
154 assert(phase_ == INIT_PHASE);
155 myPL_.set( "Max Locked", val);
156 }
157
160 if(phase_ != SOLVED_PHASE) {
161 int res = solve();
162 assert(res == 0);
163 }
164
165 return evs_;
166 }
167
170 {
171 if(phase_ != SOLVED_PHASE) {
172 int res = solve();
173 assert(res == 0);
174 }
175
176 return efs_;
177 }
178
180 virtual uint converged() const
181 {
182 if(phase_ != SOLVED_PHASE)
183 return 0;
184
185 const Anasazi::Eigensolution<ScalarT, MV>& sol = myProblem_->getSolution();
186 return sol.numVecs;
187 }
188
190 virtual uint iterations() const {
191 if(phase_ != SOLVED_PHASE)
192 return -1;
193
194 return solver_->getNumIters();
195 }
196
197 Teuchos::RCP< Anasazi::BasicEigenproblem<ScalarT, MV, OP> > getProblem()
198 {
199 return myProblem_;
200 }
201
202public:
203 virtual std::ostream& (std::ostream& os) const {
204 return os << concepts::typeOf(*this) << " in state " << (int)phase_;
205 }
206
207private:
208 Teuchos::RCP<OP_C> A_;
209 Teuchos::RCP< Anasazi::BasicEigenproblem<ScalarT, MV, OP> > myProblem_;
210
211 Teuchos::ParameterList myPL_;
212
213 Teuchos::RCP< Anasazi::SolverManager<ScalarT, MV, OP> > solver_;
214
217
218 State phase_;
219 SolverType solverType_;
220};
221
222}
void setOrthogonalization(std::string mode="SVQB")
Definition anasaziES.hh:93
virtual uint iterations() const
Returns the number of iterations.
Definition anasaziES.hh:190
virtual const concepts::Array< ScalarT > & getEV()
Returns an array with the eigen values.
Definition anasaziES.hh:159
void setConvTol(double tolerance)
Definition anasaziES.hh:129
void setMaxLocked(int val)
Definition anasaziES.hh:153
virtual uint converged() const
Returns the number of converged eigen pairs.
Definition anasaziES.hh:180
void setRelConvTol(bool rct)
Definition anasaziES.hh:138
void setConvergenceNorm(std::string norm="2")
Definition anasaziES.hh:145
virtual const concepts::Array< concepts::Vector< ScalarT > * > & getEF()
Returns an array with the eigen functions.
Definition anasaziES.hh:169
void setNumBlocks(int numBlocks)
Definition anasaziES.hh:116
std::string typeOf(const T &t)
Definition output.hh:43
Set< F > makeSet(uint n, const F &first,...)
Definition set.hh:320