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 using namespace std;
42 
43 namespace Nektar
44 {
45 namespace SolverUtils
46 {
47 
48 string DriverModifiedArnoldi::className =
49  GetDriverFactory().RegisterCreatorFunction("ModifiedArnoldi",
50  DriverModifiedArnoldi::create);
51 string DriverModifiedArnoldi::driverLookupId =
52  LibUtilities::SessionReader::RegisterEnumValue("Driver",
53  "ModifiedArnoldi",0);
54 
55 /**
56  *
57  */
58 DriverModifiedArnoldi::DriverModifiedArnoldi(
60  : DriverArnoldi(pSession)
61 {
62 }
63 
64 
65 /**
66  *
67  */
69 {
70 }
71 
72 
73 /**
74  *
75  */
77 {
79 
80  m_equ[0]->PrintSummary(out);
81 
82  // Print session parameters
83  if (m_comm->GetRank() == 0)
84  {
85  out << "\tArnoldi solver type : Modified Arnoldi" << endl;
86  }
87 
89 
90  m_equ[m_nequ - 1]->DoInitialise();
91 
92  //FwdTrans Initial conditions to be in Coefficient Space
93  m_equ[m_nequ-1] ->TransPhysToCoeff();
94 
95 }
96 
97 
98 /**
99  *
100  */
102 {
103  int i = 0;
104  int j = 0;
105  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
106  int ntot = m_nfields*nq;
107  int converged = 0;
108  NekDouble resnorm = 0.0;
109  ofstream evlout;
110  std::string evlFile = m_session->GetSessionName() + ".evl";
111 
112  if (m_comm->GetRank() == 0)
113  {
114  evlout.open(evlFile.c_str());
115  }
116 
117  // Allocate memory
122 
124  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
126  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
127  for (i = 0; i < m_kdim + 1; ++i)
128  {
129  Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
130  Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
131  }
132 
133  // Copy starting vector into second sequence element (temporary).
134  if(m_session->DefinesFunction("InitialConditions"))
135  {
136  if (m_comm->GetRank() == 0)
137  {
138  out << "\tInital vector : specified in input file " << endl;
139  }
140  m_equ[0]->SetInitialConditions(0.0,false);
141 
142  CopyFieldToArnoldiArray(Kseq[1]);
143  }
144  else
145  {
146  if (m_comm->GetRank() == 0)
147  {
148  out << "\tInital vector : random " << endl;
149  }
150 
151  NekDouble eps=0.0001;
152  Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
153  }
154 
155  // Perform one iteration to enforce boundary conditions.
156  // Set this as the initial value in the sequence.
157  EV_update(Kseq[1], Kseq[0]);
158  if (m_comm->GetRank() == 0)
159  {
160  out << "Iteration: " << 0 << endl;
161  }
162 
163  // Normalise first vector in sequence
164  alpha[0] = Blas::Ddot(ntot, &Kseq[0][0], 1, &Kseq[0][0], 1);
165  m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
166  alpha[0] = std::sqrt(alpha[0]);
167  Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
168 
169  // Fill initial krylov sequence
170  NekDouble resid0;
171  for (i = 1; !converged && i <= m_kdim; ++i)
172  {
173  // Compute next vector
174  EV_update(Kseq[i-1], Kseq[i]);
175 
176  // Normalise
177  alpha[i] = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
178  m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
179  alpha[i] = std::sqrt(alpha[i]);
180 
181  //alpha[i] = std::sqrt(alpha[i]);
182  Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
183 
184  // Copy Krylov sequence into temporary storage
185  for (int k = 0; k < i + 1; ++k)
186  {
187  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
188  }
189 
190  // Generate Hessenberg matrix and compute eigenvalues of it.
191  EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
192 
193  // Test for convergence.
194  converged = EV_test(i, i, zvec, wr, wi, resnorm,
195  std::min(i, m_nvec), evlout, resid0);
196  if ( i >= m_nvec)
197  {
198  converged = max (converged, 0);
199  }
200  else
201  {
202  converged = 0;
203  }
204 
205  if (m_comm->GetRank() == 0)
206  {
207  out << "Iteration: " << i << " (residual : " << resid0
208  << ")" <<endl;
209  }
210  }
211 
212  // Continue with full sequence
213  if (!converged)
214  {
215  for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
216  {
217  // Shift all the vectors in the sequence.
218  // First vector is removed.
219  for (int j = 1; j <= m_kdim; ++j)
220  {
221  alpha[j-1] = alpha[j];
222  Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
223  }
224 
225  // Compute next vector
226  EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
227 
228  // Compute new scale factor
229  alpha[m_kdim] = Blas::Ddot(ntot, &Kseq[m_kdim][0], 1,
230  &Kseq[m_kdim][0], 1);
231  m_comm->AllReduce(alpha[m_kdim], Nektar::LibUtilities::ReduceSum);
232  alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
233  Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
234  Kseq[m_kdim], 1);
235 
236  // Copy Krylov sequence into temporary storage
237  for (int k = 0; k < m_kdim + 1; ++k)
238  {
239  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
240  }
241 
242  // Generate Hessenberg matrix and compute eigenvalues of it
243  EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
244 
245  // Test for convergence.
246  converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm,
247  m_nvec, evlout, resid0);
248 
249  if (m_comm->GetRank() == 0)
250  {
251  out << "Iteration: " << i << " (residual : "
252  << resid0 << ")" <<endl;
253  }
254  }
255  }
256 
257  m_equ[0]->Output();
258 
259  // Evaluate and output computation time and solution accuracy.
260  // The specific format of the error output is essential for the
261  // regression tests to work.
262  // Evaluate L2 Error
263  for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
264  {
265  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
266  NekDouble vLinfError = m_equ[0]->LinfError(j);
267  if (m_comm->GetRank() == 0)
268  {
269  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
270  << ") : " << vL2Error << endl;
271  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
272  << ") : " << vLinfError << endl;
273  }
274  }
275 
276  // Process eigenvectors and write out.
277  EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
278  converged);
279 
280  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
281  "Need to implement Ritz re-evaluation of"
282  "eigenvalue. Only one half of complex "
283  "value will be correct");
284 
285  // store eigenvalues so they can be accessed from driver class
286  m_real_evl = wr;
287  m_imag_evl = wi;
288 
289  // Close the runtime info file.
290  if (m_comm->GetRank() == 0)
291  {
292  evlout.close();
293  }
294 }
295 
296 
297 /**
298  *
299  */
303 {
304  // Copy starting vector into first sequence element.
306  m_equ[0]->TransCoeffToPhys();
307 
308  m_equ[0]->DoSolve();
309 
311  {
313  fields = m_equ[0]->UpdateFields();
314 
315  //start Adjoint with latest fields of direct
316  CopyFwdToAdj();
317  m_equ[1]->TransCoeffToPhys();
318 
319  m_equ[1]->DoSolve();
320  }
321 
322  // Copy starting vector into first sequence element.
324 }
325 
326 
327 /**
328  * Computes the Ritz eigenvalues and eigenvectors of the Hessenberg matrix
329  * constructed from the Krylov sequence.
330  */
333  const int ntot,
334  const Array<OneD, NekDouble> &alpha,
335  const int kdim,
339  NekDouble &resnorm)
340 {
341  int kdimp = kdim + 1;
342  int lwork = 10*kdim;
343  int ier;
344  Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
345  Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
346  Array<OneD, NekDouble> rwork(lwork, 0.0);
347 
348  // Modified G-S orthonormalisation
349  for (int i = 0; i < kdimp; ++i)
350  {
351  NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
352  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
353  gsc = std::sqrt(gsc);
354  ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
355 
356  R[i*kdimp+i] = gsc;
357  Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
358 
359  for (int j = i + 1; j < kdimp; ++j)
360  {
361  gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
362  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
363  Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
364  R[j*kdimp+i] = gsc;
365  }
366  }
367 
368  // Compute H matrix
369  for (int i = 0; i < kdim; ++i)
370  {
371  for (int j = 0; j < kdim; ++j)
372  {
373  H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
374  - Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
375  H[j*kdim+i] /= R[j*kdimp+j];
376  }
377  }
378 
379  H[(kdim-1)*kdim+kdim] = alpha[kdim]
380  * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
381 
382  Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
383  &zvec[0], kdim, &rwork[0], lwork, ier);
384 
385  ASSERTL0(!ier, "Error with dgeev");
386 
387  resnorm = H[(kdim-1)*kdim + kdim];
388 }
389 
390 
391 /**
392  * Check for convergence of the residuals of the eigenvalues of H.
393  */
395  const int itrn,
396  const int kdim,
400  const NekDouble resnorm,
401  const int nvec,
402  ofstream &evlout,
403  NekDouble &resid0)
404 {
405  int idone = 0;
406  // NekDouble period = 0.1;
407 
408  Array<OneD, NekDouble> resid(kdim);
409  for (int i = 0; i < kdim; ++i)
410  {
411  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
412  &zvec[0] + i*kdim, 1));
413  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
414  if (wi[i] < 0.0)
415  {
416  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
417  }
418  }
419 
420  EV_sort(zvec, wr, wi, resid, kdim);
421 
422  if (resid[nvec-1] < m_evtol)
423  {
424  idone = nvec;
425  }
426 
427  if (m_comm->GetRank() == 0)
428  {
429  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
430  << resnorm << endl;
431  evlout.precision(4);
432  evlout.setf(ios::scientific, ios::floatfield);
434  {
435  evlout << " Magnitude Angle Growth "
436  << "Frequency Residual" << endl;
437  }
438  else
439  {
440  evlout << " Real Imaginary inverse real "
441  << "inverse imag Residual" << endl;
442  }
443 
444  for (int i = 0; i < kdim; i++)
445  {
446  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
447  }
448  }
449 
450  resid0 = resid[nvec-1];
451  return idone;
452 }
453 
454 
455 /**
456  * Sorts the computed eigenvalues by smallest residual first.
457  */
463  const int dim)
464 {
465  Array<OneD, NekDouble> z_tmp(dim,0.0);
466  NekDouble wr_tmp, wi_tmp, te_tmp;
467  for (int j = 1; j < dim; ++j)
468  {
469  wr_tmp = wr[j];
470  wi_tmp = wi[j];
471  te_tmp = test[j];
472  Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
473  int i = j - 1;
474  while (i >= 0 && test[i] > te_tmp)
475  {
476  wr[i+1] = wr[i];
477  wi[i+1] = wi[i];
478  test[i+1] = test[i];
479  Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
480  i--;
481  }
482  wr[i+1] = wr_tmp;
483  wi[i+1] = wi_tmp;
484  test[i+1] = te_tmp;
485  Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
486  }
487 }
488 
489 
490 /**
491  * Post-process the Ritz eigenvalues/eigenvectors of the matrix H, to compute
492  * estimations of the leading eigenvalues and eigenvectors of the original
493  * matrix.
494  */
498  const int ntot,
499  const int kdim,
500  const int nvec,
504  const int icon)
505 {
506  if (icon == 0)
507  {
508  // Not converged, write final Krylov vector
509  ASSERTL0(false, "Convergence was not achieved within the "
510  "prescribed number of iterations.");
511  }
512  else if (icon < 0)
513  {
514  // Minimum residual reached
515  ASSERTL0(false, "Minimum residual reached.");
516  }
517  else if (icon == nvec)
518  {
519  // Converged, write out eigenvectors
520  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
522  = m_equ[0]->UpdateFields();
523 
524  for (int j = 0; j < icon; ++j)
525  {
526  std::string file = m_session->GetSessionName() + "_eig_"
527  + boost::lexical_cast<std::string>(j)
528  + ".fld";
529 
530  WriteEvs(cout, j, wr[j], wi[j]);
531  WriteFld(file,Kseq[j]);
532  }
533  }
534  else
535  {
536  // Not recognised
537  ASSERTL0(false, "Unrecognised value.");
538  }
539 }
540 
541 
542 /**
543  * Compute estimates of the eigenvalues/eigenvectors of the original system.
544  */
548  const int ntot,
549  const int kdim,
550  const int nvec,
554 {
555  NekDouble wgt, norm;
556 
557  // Generate big eigenvectors
558  for (int j = 0; j < nvec; ++j)
559  {
560  Vmath::Zero(ntot, evecs[j], 1);
561  for (int i = 0; i < kdim; ++i)
562  {
563  wgt = zvec[i + j*kdim];
564  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
565  }
566  }
567 
568  // Normalise the big eigenvectors
569  for (int i = 0; i < nvec; ++i)
570  {
571  if (wi[i] == 0.0) // Real mode
572  {
573  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
574  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
575  norm = std::sqrt(norm);
576  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
577  }
578  else
579  {
580  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
581  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
582  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
583  norm = std::sqrt(norm);
584 
585  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
586  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
587 
588  i++;
589  }
590  }
591 }
592 
593 }
594 }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:70
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
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, std::ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
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)
STL namespace.
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.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
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.
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
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:194
virtual void v_InitObject(std::ostream &out=std::cout)
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
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:66
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
virtual void v_Execute(std::ostream &out=std::cout)
Virtual function for solve implementation.
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
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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