Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Static Protected Attributes | Static Private Attributes | Friends | List of all members
Nektar::SolverUtils::DriverArpack Class Reference

Base class for the development of solvers. More...

#include <DriverArpack.h>

Inheritance diagram for Nektar::SolverUtils::DriverArpack:
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className = GetDriverFactory().RegisterCreatorFunction("Arpack", DriverArpack::create)
 Name of the class. More...
 

Protected Member Functions

 DriverArpack (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
virtual ~DriverArpack ()
 Destructor. More...
 
virtual void v_InitObject (std::ostream &out=std::cout)
 Virtual function for initialisation implementation. More...
 
virtual void v_Execute (std::ostream &out=std::cout)
 Virtual function for solve implementation. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
 DriverArnoldi (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
virtual ~DriverArnoldi ()
 Destructor. More...
 
void CopyArnoldiArrayToField (Array< OneD, NekDouble > &array)
 Copy Arnoldi storage to fields. More...
 
void CopyFieldToArnoldiArray (Array< OneD, NekDouble > &array)
 Copy fields to Arnoldi storage. More...
 
void CopyFwdToAdj ()
 Copy the forward field to the adjoint system in transient growth calculations. More...
 
void WriteFld (std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
 Write coefficients to file. More...
 
void WriteFld (std::string file, Array< OneD, NekDouble > coeffs)
 
void WriteEvs (std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
 
void MaskInit ()
 init mask More...
 
void GetUnmaskFunction (std::vector< std::vector< LibUtilities::EquationSharedPtr > > &unmaskfun)
 
virtual Array< OneD, NekDoublev_GetRealEvl (void)
 
virtual Array< OneD, NekDoublev_GetImagEvl (void)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Initialises EquationSystem class members. More...
 

Protected Attributes

int m_maxn
 
int m_maxnev
 
int m_maxncv
 
- Protected Attributes inherited from Nektar::SolverUtils::DriverArnoldi
int m_kdim
 
int m_nvec
 Dimension of Krylov subspace. More...
 
int m_nits
 Number of vectors to test. More...
 
NekDouble m_evtol
 Maxmum number of iterations. More...
 
NekDouble m_period
 Tolerance of iteratiosn. More...
 
bool m_timeSteppingAlgorithm
 Period of time stepping algorithm. More...
 
int m_infosteps
 underlying operator is time stepping More...
 
int m_nfields
 interval to dump information if required. More...
 
NekDouble m_realShift
 
NekDouble m_imagShift
 
int m_negatedOp
 
Array< OneD, NekDoublem_real_evl
 Operator in solve call is negated. More...
 
Array< OneD, NekDoublem_imag_evl
 
bool m_useMask
 
Array< OneD, NekDoublem_maskCoeffs
 
Array< OneD, NekDoublem_maskPhys
 
- Protected Attributes inherited from Nektar::SolverUtils::Driver
LibUtilities::CommSharedPtr m_comm
 Communication object. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session reader object. More...
 
LibUtilities::SessionReaderSharedPtr session_LinNS
 I the Coupling between SFD and arnoldi. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 MeshGraph object. More...
 
Array< OneD, EquationSystemSharedPtrm_equ
 Equation system to solve. More...
 
int m_nequ
 number of equations More...
 
enum EvolutionOperatorType m_EvolutionOperator
 Evolution Operator. More...
 

Static Protected Attributes

static std::string arpackProblemTypeLookupIds []
 
static std::string arpackProblemTypeDefault = LibUtilities::SessionReader::RegisterDefaultSolverInfo("ArpackProblemType","LargestMag")
 
static std::string driverLookupId = LibUtilities::SessionReader::RegisterEnumValue("Driver","Arpack",0)
 
- Static Protected Attributes inherited from Nektar::SolverUtils::Driver
static std::string evolutionOperatorLookupIds []
 
static std::string evolutionOperatorDef
 
static std::string driverDefault
 

Static Private Attributes

static std::string ArpackProblemTypeTrans []
 

Friends

class MemoryManager< DriverArpack >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
SOLVER_UTILS_EXPORT void ArnoldiSummary (std::ostream &out)
 
SOLVER_UTILS_EXPORT const Array< OneD, const NekDouble > & GetMaskCoeff () const
 
SOLVER_UTILS_EXPORT const Array< OneD, const NekDouble > & GetMaskPhys () const
 
- Public Member Functions inherited from Nektar::SolverUtils::Driver
virtual ~Driver ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (std::ostream &out=std::cout)
 Initialise Object. More...
 
SOLVER_UTILS_EXPORT void Execute (std::ostream &out=std::cout)
 Execute driver. More...
 
SOLVER_UTILS_EXPORT Array< OneD, EquationSystemSharedPtrGetEqu ()
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetRealEvl (void)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetImagEvl (void)
 

Detailed Description

Base class for the development of solvers.

Definition at line 48 of file DriverArpack.h.

Constructor & Destructor Documentation

◆ DriverArpack()

Nektar::SolverUtils::DriverArpack::DriverArpack ( const LibUtilities::SessionReaderSharedPtr  pSession,
const SpatialDomains::MeshGraphSharedPtr  pGraph 
)
protected

Constructor.

Definition at line 64 of file DriverArpack.cpp.

66  : DriverArnoldi(pSession, pGraph)
67 {
68 }
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.

◆ ~DriverArpack()

Nektar::SolverUtils::DriverArpack::~DriverArpack ( )
protectedvirtual

Destructor.

Definition at line 74 of file DriverArpack.cpp.

75 {
76 }

Member Function Documentation

◆ create()

static DriverSharedPtr Nektar::SolverUtils::DriverArpack::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 54 of file DriverArpack.h.

57  {
59  pSession, pGraph);
60  p->InitObject();
61  return p;
62  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:51

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_Execute()

void Nektar::SolverUtils::DriverArpack::v_Execute ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Definition at line 119 of file DriverArpack.cpp.

121 {
122  Array<OneD, NekDouble> tmpworkd;
123 
124  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs(); // Number of points in the mesh
125  int n = m_nfields*nq; // Number of points in eigenvalue calculation
126  int lworkl = 3*m_kdim*(m_kdim+2); // Size of work array
127  int ido ; //REVERSE COMMUNICATION parameter. At the first call must be initialised at 0
128  int info; // do not set initial vector (info=0 random initial vector, info=1 read initial vector from session file)
129 
130  int iparam[11];
131  int ipntr[14];
132 
133  Array<OneD, int> ritzSelect;
134  Array<OneD, NekDouble> dr;
135  Array<OneD, NekDouble> di;
136  Array<OneD, NekDouble> workev;
137  Array<OneD, NekDouble> z;
138  NekDouble sigmar, sigmai;
139 
140  Array<OneD, NekDouble> resid(n);
141  Array<OneD, NekDouble> v(n*m_kdim);
142  Array<OneD, NekDouble> workl(lworkl, 0.0);
143  Array<OneD, NekDouble> workd(3*n, 0.0);
144 
145  ASSERTL0(n <= m_maxn, "N is greater than MAXN");
146 
147  if(m_session->DefinesFunction("InitialConditions"))
148  {
149  out << "\tInital vector : input file " << endl;
150  info = 1;
152  }
153  else
154  {
155  out << "\tInital vector : random " << endl;
156  info = 0;
157  }
158 
159  char B;
160 
161  iparam[0] = 1; // strategy for shift-invert
162  iparam[1] = 0; // (deprecated)
163  iparam[2] = m_nits; // maximum number of iterations allowed/taken
164  iparam[3] = 1; // blocksize to be used for recurrence
165  iparam[4] = 0; // number of converged ritz eigenvalues
166  iparam[5] = 0; // (deprecated)
167 
168  // Use generalized B matrix for coupled solver.
170  {
171  iparam[6] = 1; // computation mode 1=> matrix-vector prod
172  B = 'I';
173  }
174  else {
175  iparam[6] = 3; // computation mode 1=> matrix-vector prod
176  B = 'G';
177  }
178 #if 0
179  if((fabs(m_realShift) > NekConstants::kNekZeroTol)|| // use shift if m_realShift > 1e-12
181  {
182  iparam[6] = 3; // This was 3 need to know what to set it to
183  B = 'G';
184  }
185  else
186  {
187  iparam[6] = 1; // computation mode 1=> matrix-vector prod
188  B = 'I';
189  }
190 #endif
191  iparam[7] = 0; // (for shift-invert)
192  iparam[8] = 0; // number of MV operations
193  iparam[9] = 0; // number of BV operations
194  iparam[10]= 0; // number of reorthogonalisation steps
195 
196  int cycle = 0;
197  const char* problem = ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")].c_str();
198 
199  std::string name = m_session->GetSessionName() + ".evl";
200  ofstream pFile(name.c_str());
201 
202  ido = 0; //At the first call must be initialisedat 0
203 
204  while(ido != 99)//ido==-1 || ido==1 || ido==0)
205  {
206  //Routine for eigenvalue evaluation for non-symmetric operators
207  Arpack::Dnaupd( ido, &B, // B='I' for std eval problem
208  n, problem, m_nvec,
209  m_evtol, &resid[0], m_kdim,
210  &v[0], n, iparam, ipntr, &workd[0],
211  &workl[0], lworkl, info);
212 
213  //Plotting of real and imaginary part of the
214  //eigenvalues from workl
215  out << "\rIteration " << cycle << ", output: " << info << ", ido=" << ido << " " << std::flush;
216 
217  if(!((cycle-1)%m_kdim)&&(cycle> m_kdim)&&(ido!=2))
218  {
219  pFile << "Krylov spectrum at iteration: " << cycle << endl;
220 
222  {
223  pFile << "EV Magnitude Angle Growth Frequency Residual" << endl;
224  }
225  else
226  {
227  pFile << "EV Real Imaginary inverse real inverse imag Residual" << endl;
228  }
229 
230  out << endl;
231  for(int k = 0; k < m_kdim; ++k)
232  {
233  // write m_kdim eigs to screen
234  WriteEvs(pFile,k, workl[ipntr[5]-1+k],workl[ipntr[6]-1+k]);
235  }
236  }
237 
238  if (ido == 99) break;
239 
240  switch(ido)
241  {
242  case -1:
243  case 1: // Note that ido=1 we are using input x
244  // (workd[inptr[0]-1]) rather than Mx as
245  // recommended in manual since it is not
246  // possible to impose forcing directly.
247  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0]-1));
248 
249  m_equ[0]->TransCoeffToPhys();
250 
251  m_equ[0]->DoSolve();
253  {
254  //start Adjoint with latest fields of direct
255  CopyFwdToAdj();
256 
257  m_equ[1]->TransCoeffToPhys();
258  m_equ[1]->DoSolve();
259  }
260 
261  if(!(cycle%m_infosteps))
262  {
263  out << endl;
264  m_equ[0]->Output();
265  }
266 
267  // operated fields are copied into workd[inptr[1]-1]
268  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1]-1));
269 
270  cycle++;
271  break;
272  case 2: // provide y = M x (bwd trans and iproduct);
273  {
274  //workd[inptr[0]-1] copied into operator fields
275  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0]-1));
276 
277  m_equ[0]->TransCoeffToPhys();
278 
279  Array<OneD, MultiRegions::ExpListSharedPtr> fields = m_equ[0]->UpdateFields();
280  for (int i = 0; i < fields.size(); ++i)
281  {
282  fields[i]->IProductWRTBase(fields[i]->GetPhys(),
283  fields[i]->UpdateCoeffs());
284  }
285 
286  // operated fields are copied into workd[inptr[1]-1]
287  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1]-1));
288  break;
289  }
290  default:
291  ASSERTL0(false, "Unexpected reverse communication request.");
292  }
293 
294  }
295 
296  out<< endl << "Converged in " << iparam[8] << " iterations" << endl;
297 
298  ASSERTL0(info >= 0," Error with Dnaupd");
299 
300  ritzSelect = Array<OneD, int> (m_kdim,0);
301  dr = Array<OneD, NekDouble> (m_nvec+1,0.0);
302  di = Array<OneD, NekDouble> (m_nvec+1,0.0);
303  workev = Array<OneD, NekDouble> (3*m_kdim);
304  z = Array<OneD, NekDouble> (n*(m_nvec+1));
305 
306  if(m_negatedOp)
307  {
308  sigmar = -m_realShift;
309  }
310  else
311  {
312  sigmar = m_realShift;
313  }
314 
315  // Do not pass imaginary shift to Arpack since we have not
316  // used a Fortran complex number format and so processing
317  // is mucked up. Need to do some processing afterwards.
318  sigmai = 0;
319 
320  //Setting 'A', Ritz vectors are computed. 'S' for Shur vectors
321  Arpack::Dneupd(1, "A", ritzSelect.get(), dr.get(), di.get(),
322  z.get(), n, sigmar, sigmai, workev.get(), &B,
323  n, problem, m_nvec, m_evtol, resid.get(), m_kdim,
324  v.get(), n, iparam, ipntr, workd.get(),
325  workl.get(),lworkl,info);
326 
327  ASSERTL0(info == 0, " Error with Dneupd");
328 
329  int nconv=iparam[4];
330 
331  // Subtract off complex shift if it exists
332  if(m_negatedOp)
333  {
334  Vmath::Sadd(nconv,m_imagShift,di,1,di,1);
335  }
336  else
337  {
338  Vmath::Sadd(nconv,-m_imagShift,di,1,di,1);
339  }
340 
341  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
342  "Need to implement Ritz re-evaluation of"
343  "eigenvalue. Only one half of complex "
344  "value will be correct");
345 
346 
347  Array<OneD, MultiRegions::ExpListSharedPtr> fields = m_equ[0]->UpdateFields();
348 
349  out << "Converged Eigenvalues: " << nconv << endl;
350  pFile << "Converged Eigenvalues: " << nconv << endl;
351 
353  {
354  out << " Magnitude Angle Growth Frequency" << endl;
355  pFile << " Magnitude Angle Growth Frequency" << endl;
356  for(int i= 0; i< nconv; ++i)
357  {
358  WriteEvs(out,i,dr[i],di[i]);
359  WriteEvs(pFile,i,dr[i],di[i]);
360 
361  std::string file = m_session->GetSessionName() + "_eig_"
362  + boost::lexical_cast<std::string>(i)
363  + ".fld";
364  WriteFld(file,z + i*n);
365  }
366  }
367  else
368  {
369  out << " Real Imaginary " << endl;
370  pFile << " Real Imaginary " << endl;
371  for(int i= 0; i< nconv; ++i)
372  {
373  WriteEvs(out,i,dr[i],di[i],
375  WriteEvs(pFile,i,dr[i],di[i],
377 
378  std::string file = m_session->GetSessionName() + "_eig_"
379  + boost::lexical_cast<std::string>(i)
380  + ".fld";
381  WriteFld(file,z + i*n);
382  }
383  }
384 
385  m_real_evl = dr;
386  m_imag_evl = di;
387 
388  pFile.close();
389 
390  for(int j = 0; j < m_equ[0]->GetNvariables(); ++j)
391  {
392  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
393  NekDouble vLinfError = m_equ[0]->LinfError(j);
394  if (m_comm->GetRank() == 0)
395  {
396  out << "L 2 error (variable " << m_equ[0]->GetVariable(j) << ") : " << vL2Error << endl;
397  out << "L inf error (variable " << m_equ[0]->GetVariable(j) << ") : " << vLinfError << endl;
398  }
399  }
400 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:223
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:65
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:59
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:63
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:60
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:73
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:61
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:67
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:72
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
static std::string ArpackProblemTypeTrans[]
Definition: DriverArpack.h:92
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:89
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:104
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98
static void Dnaupd(int &ido, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Top level reverse communication interface to solve real double-precision non-symmetric problems.
Definition: Arpack.hpp:134
static void Dneupd(const int &rvec, const char *howmny, const int *select, double *dr, double *di, double *z, const int &ldz, const double &sigmar, const double &sigmai, double *workev, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Post-processing routine to computed eigenvector of computed eigenvalues in Dnaupd.
Definition: Arpack.hpp:148
static const NekDouble kNekUnsetDouble
static const NekDouble kNekZeroTol
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
Definition: Vmath.cpp:341

References ArpackProblemTypeTrans, ASSERTL0, Nektar::SolverUtils::DriverArnoldi::CopyArnoldiArrayToField(), Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Nektar::SolverUtils::DriverArnoldi::CopyFwdToAdj(), Arpack::Dnaupd(), Arpack::Dneupd(), Nektar::SolverUtils::eTransientGrowth, Nektar::NekConstants::kNekUnsetDouble, Nektar::NekConstants::kNekZeroTol, Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_EvolutionOperator, Nektar::SolverUtils::DriverArnoldi::m_evtol, Nektar::SolverUtils::DriverArnoldi::m_imag_evl, Nektar::SolverUtils::DriverArnoldi::m_imagShift, Nektar::SolverUtils::DriverArnoldi::m_infosteps, Nektar::SolverUtils::DriverArnoldi::m_kdim, m_maxn, Nektar::SolverUtils::DriverArnoldi::m_negatedOp, Nektar::SolverUtils::DriverArnoldi::m_nfields, Nektar::SolverUtils::DriverArnoldi::m_nits, Nektar::SolverUtils::DriverArnoldi::m_nvec, Nektar::SolverUtils::DriverArnoldi::m_real_evl, Nektar::SolverUtils::DriverArnoldi::m_realShift, Nektar::SolverUtils::Driver::m_session, Nektar::SolverUtils::DriverArnoldi::m_timeSteppingAlgorithm, CellMLToNektar.pycml::name, Vmath::Sadd(), WARNINGL0, Nektar::SolverUtils::DriverArnoldi::WriteEvs(), and Nektar::SolverUtils::DriverArnoldi::WriteFld().

◆ v_InitObject()

void Nektar::SolverUtils::DriverArpack::v_InitObject ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Definition at line 81 of file DriverArpack.cpp.

82 {
84 
85  //Initialisation of Arnoldi parameters
86  m_maxn = 1000000; // Maximum size of the problem
87  m_maxnev = 200; // maximum number of eigenvalues requested
88  m_maxncv = 500; // Largest number of basis vector used in Implicitly Restarted Arnoldi
89 
90  // Error alerts
91  ASSERTL0(m_nvec <= m_maxnev,"NEV is greater than MAXNEV");
92  ASSERTL0(m_kdim <= m_maxncv,"NEV is greater than MAXNEV");
93  ASSERTL0(2 <= m_kdim-m_nvec,"NCV-NEV is less than 2");
94 
95  ASSERTL0(m_comm->GetSize()==1,"..ARPACK is not currently set-up for parallel execution...\n");
96  m_equ[0]->PrintSummary(out);
97 
98  ASSERTL0(m_comm->GetSize() == 1,
99  "ARPACK Arnoldi solver does not support execution in parallel.");
100 
101  // Print session parameters
102  out << "\tArnoldi solver type : Arpack" << endl;
103 
104  out << "\tArpack problem type : ";
105  out << ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")] << endl;
107 
108  for( int i = 0; i < m_nequ; ++i)
109  {
110  m_equ[i]->DoInitialise();
111  }
112 
113  // FwdTrans Initial conditions to be in Coefficient Space
114  m_equ[m_nequ-1] ->TransPhysToCoeff();
115 
116 }
virtual void v_InitObject(std::ostream &out=std::cout)
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
int m_nequ
number of equations
Definition: Driver.h:101

References Nektar::SolverUtils::DriverArnoldi::ArnoldiSummary(), ArpackProblemTypeTrans, ASSERTL0, Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::DriverArnoldi::m_kdim, m_maxn, m_maxncv, m_maxnev, Nektar::SolverUtils::Driver::m_nequ, Nektar::SolverUtils::DriverArnoldi::m_nvec, Nektar::SolverUtils::Driver::m_session, and Nektar::SolverUtils::DriverArnoldi::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< DriverArpack >

friend class MemoryManager< DriverArpack >
friend

Definition at line 1 of file DriverArpack.h.

Member Data Documentation

◆ arpackProblemTypeDefault

std::string Nektar::SolverUtils::DriverArpack::arpackProblemTypeDefault = LibUtilities::SessionReader::RegisterDefaultSolverInfo("ArpackProblemType","LargestMag")
staticprotected

Definition at line 88 of file DriverArpack.h.

◆ arpackProblemTypeLookupIds

std::string Nektar::SolverUtils::DriverArpack::arpackProblemTypeLookupIds
staticprotected
Initial value:
= {
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestReal" ,0),
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestReal" ,1),
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestImag" ,2),
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestImag" ,3),
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestMag" ,4),
LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestMag" ,5),
}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 87 of file DriverArpack.h.

◆ ArpackProblemTypeTrans

std::string Nektar::SolverUtils::DriverArpack::ArpackProblemTypeTrans
staticprivate
Initial value:
=
{ "LR", "SR", "LI", "SI", "LM", "SM" }

Definition at line 92 of file DriverArpack.h.

Referenced by v_Execute(), and v_InitObject().

◆ className

std::string Nektar::SolverUtils::DriverArpack::className = GetDriverFactory().RegisterCreatorFunction("Arpack", DriverArpack::create)
static

Name of the class.

Definition at line 65 of file DriverArpack.h.

◆ driverLookupId

std::string Nektar::SolverUtils::DriverArpack::driverLookupId = LibUtilities::SessionReader::RegisterEnumValue("Driver","Arpack",0)
staticprotected

Definition at line 89 of file DriverArpack.h.

◆ m_maxn

int Nektar::SolverUtils::DriverArpack::m_maxn
protected

Definition at line 70 of file DriverArpack.h.

Referenced by v_Execute(), and v_InitObject().

◆ m_maxncv

int Nektar::SolverUtils::DriverArpack::m_maxncv
protected

Definition at line 72 of file DriverArpack.h.

Referenced by v_InitObject().

◆ m_maxnev

int Nektar::SolverUtils::DriverArpack::m_maxnev
protected

Definition at line 71 of file DriverArpack.h.

Referenced by v_InitObject().