Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DriverModifiedArnoldi.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File DriverModifiedArnoldi.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Driver for eigenvalue analysis using the modified Arnoldi
33 // method.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
37 #include <iomanip>
38 
40 
41 namespace Nektar
42 {
43 namespace SolverUtils
44 {
45 
47  GetDriverFactory().RegisterCreatorFunction("ModifiedArnoldi",
51  "ModifiedArnoldi",0);
52 
53 /**
54  *
55  */
58  : DriverArnoldi(pSession)
59 {
60 }
61 
62 
63 /**
64  *
65  */
67 {
68 }
69 
70 
71 /**
72  *
73  */
75 {
77 
78  m_equ[0]->PrintSummary(out);
79 
80  // Print session parameters
81  if (m_comm->GetRank() == 0)
82  {
83  out << "\tArnoldi solver type : Modified Arnoldi" << endl;
84  }
85 
87 
88  m_equ[m_nequ - 1]->DoInitialise();
89 
90  //FwdTrans Initial conditions to be in Coefficient Space
91  m_equ[m_nequ-1] ->TransPhysToCoeff();
92 
93 }
94 
95 
96 /**
97  *
98  */
100 {
101  int i = 0;
102  int j = 0;
103  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
104  int ntot = m_nfields*nq;
105  int converged = 0;
106  NekDouble resnorm = 0.0;
107  ofstream evlout;
108  std::string evlFile = m_session->GetSessionName() + ".evl";
109 
110  if (m_comm->GetRank() == 0)
111  {
112  evlout.open(evlFile.c_str());
113  }
114 
115  // Allocate memory
120 
122  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
124  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
125  for (i = 0; i < m_kdim + 1; ++i)
126  {
127  Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
128  Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
129  }
130 
131  // Copy starting vector into second sequence element (temporary).
132  if(m_session->DefinesFunction("InitialConditions"))
133  {
134  if (m_comm->GetRank() == 0)
135  {
136  out << "\tInital vector : specified in input file " << endl;
137  }
138  m_equ[0]->SetInitialConditions(0.0,false);
139 
140  CopyFieldToArnoldiArray(Kseq[1]);
141  }
142  else
143  {
144  if (m_comm->GetRank() == 0)
145  {
146  out << "\tInital vector : random " << endl;
147  }
148 
149  NekDouble eps=0.0001;
150  Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
151  }
152 
153  // Perform one iteration to enforce boundary conditions.
154  // Set this as the initial value in the sequence.
155  EV_update(Kseq[1], Kseq[0]);
156  if (m_comm->GetRank() == 0)
157  {
158  out << "Iteration: " << 0 << endl;
159  }
160 
161  // Normalise first vector in sequence
162  alpha[0] = Blas::Ddot(ntot, &Kseq[0][0], 1, &Kseq[0][0], 1);
163  m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
164  alpha[0] = std::sqrt(alpha[0]);
165  Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
166 
167  // Fill initial krylov sequence
168  NekDouble resid0;
169  for (i = 1; !converged && i <= m_kdim; ++i)
170  {
171  // Compute next vector
172  EV_update(Kseq[i-1], Kseq[i]);
173 
174  // Normalise
175  alpha[i] = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
176  m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
177  alpha[i] = std::sqrt(alpha[i]);
178 
179  //alpha[i] = std::sqrt(alpha[i]);
180  Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
181 
182  // Copy Krylov sequence into temporary storage
183  for (int k = 0; k < i + 1; ++k)
184  {
185  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
186  }
187 
188  // Generate Hessenberg matrix and compute eigenvalues of it.
189  EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
190 
191  // Test for convergence.
192  converged = EV_test(i, i, zvec, wr, wi, resnorm,
193  std::min(i, m_nvec), evlout, resid0);
194  converged = max (converged, 0);
195 
196  if (m_comm->GetRank() == 0)
197  {
198  out << "Iteration: " << i << " (residual : " << resid0
199  << ")" <<endl;
200  }
201  }
202 
203  // Continue with full sequence
204  if (!converged)
205  {
206  for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
207  {
208  // Shift all the vectors in the sequence.
209  // First vector is removed.
210  for (int j = 1; j <= m_kdim; ++j)
211  {
212  alpha[j-1] = alpha[j];
213  Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
214  }
215 
216  // Compute next vector
217  EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
218 
219  // Compute new scale factor
220  alpha[m_kdim] = Blas::Ddot(ntot, &Kseq[m_kdim][0], 1,
221  &Kseq[m_kdim][0], 1);
222  m_comm->AllReduce(alpha[m_kdim], Nektar::LibUtilities::ReduceSum);
223  alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
224  Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
225  Kseq[m_kdim], 1);
226 
227  // Copy Krylov sequence into temporary storage
228  for (int k = 0; k < m_kdim + 1; ++k)
229  {
230  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
231  }
232 
233  // Generate Hessenberg matrix and compute eigenvalues of it
234  EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
235 
236  // Test for convergence.
237  converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm,
238  m_nvec, evlout, resid0);
239 
240  if (m_comm->GetRank() == 0)
241  {
242  out << "Iteration: " << i << " (residual : "
243  << resid0 << ")" <<endl;
244  }
245  }
246  }
247 
248  m_equ[0]->Output();
249 
250  // Evaluate and output computation time and solution accuracy.
251  // The specific format of the error output is essential for the
252  // regression tests to work.
253  // Evaluate L2 Error
254  for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
255  {
256  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
257  NekDouble vLinfError = m_equ[0]->LinfError(j);
258  if (m_comm->GetRank() == 0)
259  {
260  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
261  << ") : " << vL2Error << endl;
262  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
263  << ") : " << vLinfError << endl;
264  }
265  }
266 
267  // Process eigenvectors and write out.
268  EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
269  converged);
270 
271  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
272  "Need to implement Ritz re-evaluation of"
273  "eigenvalue. Only one half of complex "
274  "value will be correct");
275 
276  // store eigenvalues so they can be accessed from driver class
277  m_real_evl = wr;
278  m_imag_evl = wi;
279 
280  // Close the runtime info file.
281  if (m_comm->GetRank() == 0)
282  {
283  evlout.close();
284  }
285 }
286 
287 
288 /**
289  *
290  */
294 {
295  // Copy starting vector into first sequence element.
297  m_equ[0]->TransCoeffToPhys();
298 
299  m_equ[0]->DoSolve();
300 
302  {
304  fields = m_equ[0]->UpdateFields();
305 
306  //start Adjoint with latest fields of direct
307  CopyFwdToAdj();
308  m_equ[1]->TransCoeffToPhys();
309 
310  m_equ[1]->DoSolve();
311  }
312 
313  // Copy starting vector into first sequence element.
315 }
316 
317 
318 /**
319  * Computes the Ritz eigenvalues and eigenvectors of the Hessenberg matrix
320  * constructed from the Krylov sequence.
321  */
324  const int ntot,
325  const Array<OneD, NekDouble> &alpha,
326  const int kdim,
330  NekDouble &resnorm)
331 {
332  int kdimp = kdim + 1;
333  int lwork = 10*kdim;
334  int ier;
335  Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
336  Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
337  Array<OneD, NekDouble> rwork(lwork, 0.0);
338 
339  // Modified G-S orthonormalisation
340  for (int i = 0; i < kdimp; ++i)
341  {
342  NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
343  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
344  gsc = std::sqrt(gsc);
345  ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
346 
347  R[i*kdimp+i] = gsc;
348  Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
349 
350  for (int j = i + 1; j < kdimp; ++j)
351  {
352  gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
353  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
354  Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
355  R[j*kdimp+i] = gsc;
356  }
357  }
358 
359  // Compute H matrix
360  for (int i = 0; i < kdim; ++i)
361  {
362  for (int j = 0; j < kdim; ++j)
363  {
364  H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
365  - Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
366  H[j*kdim+i] /= R[j*kdimp+j];
367  }
368  }
369 
370  H[(kdim-1)*kdim+kdim] = alpha[kdim]
371  * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
372 
373  Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
374  &zvec[0], kdim, &rwork[0], lwork, ier);
375 
376  ASSERTL0(!ier, "Error with dgeev");
377 
378  resnorm = H[(kdim-1)*kdim + kdim];
379 }
380 
381 
382 /**
383  * Check for convergence of the residuals of the eigenvalues of H.
384  */
386  const int itrn,
387  const int kdim,
391  const NekDouble resnorm,
392  const int nvec,
393  ofstream &evlout,
394  NekDouble &resid0)
395 {
396  int idone = 0;
397  // NekDouble period = 0.1;
398 
399  Array<OneD, NekDouble> resid(kdim);
400  for (int i = 0; i < kdim; ++i)
401  {
402  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
403  &zvec[0] + i*kdim, 1));
404  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
405  if (wi[i] < 0.0)
406  {
407  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
408  }
409  }
410 
411  EV_sort(zvec, wr, wi, resid, kdim);
412 
413  if (resid[nvec-1] < m_evtol)
414  {
415  idone = nvec;
416  }
417 
418  if (m_comm->GetRank() == 0)
419  {
420  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
421  << resnorm << endl;
422  evlout.precision(4);
423  evlout.setf(ios::scientific, ios::floatfield);
425  {
426  evlout << " Magnitude Angle Growth "
427  << "Frequency Residual" << endl;
428  }
429  else
430  {
431  evlout << " Real Imaginary inverse real "
432  << "inverse imag Residual" << endl;
433  }
434 
435  for (int i = 0; i < kdim; i++)
436  {
437  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
438  }
439  }
440 
441  resid0 = resid[nvec-1];
442  return idone;
443 }
444 
445 
446 /**
447  * Sorts the computed eigenvalues by smallest residual first.
448  */
454  const int dim)
455 {
456  Array<OneD, NekDouble> z_tmp(dim,0.0);
457  NekDouble wr_tmp, wi_tmp, te_tmp;
458  for (int j = 1; j < dim; ++j)
459  {
460  wr_tmp = wr[j];
461  wi_tmp = wi[j];
462  te_tmp = test[j];
463  Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
464  int i = j - 1;
465  while (i >= 0 && test[i] > te_tmp)
466  {
467  wr[i+1] = wr[i];
468  wi[i+1] = wi[i];
469  test[i+1] = test[i];
470  Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
471  i--;
472  }
473  wr[i+1] = wr_tmp;
474  wi[i+1] = wi_tmp;
475  test[i+1] = te_tmp;
476  Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
477  }
478 }
479 
480 
481 /**
482  * Post-process the Ritz eigenvalues/eigenvectors of the matrix H, to compute
483  * estimations of the leading eigenvalues and eigenvectors of the original
484  * matrix.
485  */
489  const int ntot,
490  const int kdim,
491  const int nvec,
495  const int icon)
496 {
497  if (icon == 0)
498  {
499  // Not converged, write final Krylov vector
500  ASSERTL0(false, "Convergence was not achieved within the "
501  "prescribed number of iterations.");
502  }
503  else if (icon < 0)
504  {
505  // Minimum residual reached
506  ASSERTL0(false, "Minimum residual reached.");
507  }
508  else if (icon == nvec)
509  {
510  // Converged, write out eigenvectors
511  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
513  = m_equ[0]->UpdateFields();
514 
515  for (int j = 0; j < icon; ++j)
516  {
517  std::string file = m_session->GetSessionName() + "_eig_"
518  + boost::lexical_cast<std::string>(j)
519  + ".fld";
520 
521  WriteEvs(cout, j, wr[j], wi[j]);
522  WriteFld(file,Kseq[j]);
523  }
524  }
525  else
526  {
527  // Not recognised
528  ASSERTL0(false, "Unrecognised value.");
529  }
530 }
531 
532 
533 /**
534  * Compute estimates of the eigenvalues/eigenvectors of the original system.
535  */
539  const int ntot,
540  const int kdim,
541  const int nvec,
545 {
546  NekDouble wgt, norm;
547 
548  // Generate big eigenvectors
549  for (int j = 0; j < nvec; ++j)
550  {
551  Vmath::Zero(ntot, evecs[j], 1);
552  for (int i = 0; i < kdim; ++i)
553  {
554  wgt = zvec[i + j*kdim];
555  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
556  }
557  }
558 
559  // Normalise the big eigenvectors
560  for (int i = 0; i < nvec; ++i)
561  {
562  if (wi[i] == 0.0) // Real mode
563  {
564  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
565  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
566  norm = std::sqrt(norm);
567  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
568  }
569  else
570  {
571  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
572  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
573  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
574  norm = std::sqrt(norm);
575 
576  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
577  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
578 
579  i++;
580  }
581  }
582 }
583 
584 }
585 }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:70
virtual void v_InitObject(ostream &out=cout)
Virtual function for initialisation implementation.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:58
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
static std::string className
Name of the class.
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:471
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
void EV_small(Array< OneD, Array< OneD, NekDouble > > &Kseq, const int ntot, const Array< OneD, NekDouble > &alpha, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, NekDouble &resnorm)
Generates the upper Hessenberg matrix H and computes its eigenvalues.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
DriverModifiedArnoldi(const LibUtilities::SessionReaderSharedPtr pSession)
Constructor.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
void EV_sort(Array< OneD, NekDouble > &evec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, Array< OneD, NekDouble > &test, const int dim)
Sorts a sequence of eigenvectors/eigenvalues by magnitude.
void EV_update(Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
Generates a new vector in the sequence by applying the linear operator.
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
virtual void v_Execute(ostream &out=cout)
Virtual function for solve implementation.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:56
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:69
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:167
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
double NekDouble
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:434
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:64
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:60
void WriteEvs(ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:900
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
void EV_post(Array< OneD, Array< OneD, NekDouble > > &Tseq, Array< OneD, Array< OneD, NekDouble > > &Kseq, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const int icon)
void EV_big(Array< OneD, Array< OneD, NekDouble > > &bvecs, Array< OneD, Array< OneD, NekDouble > > &evecs, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi)
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:138
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
Base class for the development of solvers.
Definition: DriverArnoldi.h:47
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
int EV_test(const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, const int nvec, ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
virtual void v_InitObject(ostream &out=cout)
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
int m_nequ
number of equations
Definition: Driver.h:97
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:57