96 "..ARPACK is not currently set-up for parallel execution...\n");
98 m_equ[0]->PrintSummary(out);
101 out <<
"\tArnoldi solver type : Arpack" << endl;
102 out <<
"\tArpack problem type : ";
104 "ArpackProblemType")]
109 for (
int i = 0; i <
m_nequ; ++i)
111 m_equ[i]->DoInitialise();
153 if (
m_session->DefinesFunction(
"InitialConditions"))
155 out <<
"\tInital vector : input file " << endl;
161 out <<
"\tInital vector : random " << endl;
204 const char *problem =
206 "ArpackProblemType")]
210 ofstream pFile(
name.c_str());
219 iparam, ipntr, &workd[0], &workl[0], lworkl, info);
223 out <<
"\rIteration " << cycle <<
", output: " << info
224 <<
", ido=" << ido <<
" " << std::flush;
226 if (!((cycle - 1) %
m_kdim) && (cycle >
m_kdim) && (ido != 2))
228 pFile <<
"Krylov spectrum at iteration: " << cycle << endl;
232 pFile <<
"EV Magnitude Angle Growth Frequency "
238 pFile <<
"EV Real Imaginary inverse real inverse "
244 for (
int k = 0; k <
m_kdim; ++k)
247 WriteEvs(pFile, k, workl[ipntr[5] - 1 + k],
248 workl[ipntr[6] - 1 + k]);
266 m_equ[0]->TransCoeffToPhys();
274 m_equ[1]->TransCoeffToPhys();
294 m_equ[0]->TransCoeffToPhys();
297 m_equ[0]->UpdateFields();
298 for (
int i = 0; i < fields.size(); ++i)
300 fields[i]->IProductWRTBase(fields[i]->GetPhys(),
301 fields[i]->UpdateCoeffs());
309 ASSERTL0(
false,
"Unexpected reverse communication request.");
313 out << endl <<
"Converged in " << iparam[8] <<
" iterations" << endl;
315 ASSERTL0(info >= 0,
" Error with Dnaupd");
338 Arpack::Dneupd(1,
"A", ritzSelect.get(), dr.get(), di.get(),
z.get(), n,
339 sigmar, sigmai, workev.get(), &B, n, problem,
m_nvec,
341 workd.get(), workl.get(), lworkl, info);
343 ASSERTL0(info == 0,
" Error with Dneupd");
345 int nconv = iparam[4];
358 "Need to implement Ritz re-evaluation of"
359 "eigenvalue. Only one half of complex "
360 "value will be correct");
363 m_equ[0]->UpdateFields();
365 out <<
"Converged Eigenvalues: " << nconv << endl;
366 pFile <<
"Converged Eigenvalues: " << nconv << endl;
370 out <<
" Magnitude Angle Growth Frequency" << endl;
371 pFile <<
" Magnitude Angle Growth Frequency"
373 for (
int i = 0; i < nconv; ++i)
378 std::string file =
m_session->GetSessionName() +
"_eig_" +
379 std::to_string(i) +
".fld";
385 out <<
" Real Imaginary " << endl;
386 pFile <<
" Real Imaginary " << endl;
387 for (
int i = 0; i < nconv; ++i)
394 std::string file =
m_session->GetSessionName() +
"_eig_" +
395 std::to_string(i) +
".fld";
405 for (
int j = 0; j <
m_equ[0]->GetNvariables(); ++j)
409 if (
m_comm->GetRank() == 0)
411 out <<
"L 2 error (variable " <<
m_equ[0]->GetVariable(j)
412 <<
") : " << vL2Error << endl;
413 out <<
"L inf error (variable " <<
m_equ[0]->GetVariable(j)
414 <<
") : " << vLinfError << endl;
#define ASSERTL0(condition, msg)
#define WARNINGL0(condition, msg)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
Base class for the development of solvers.
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
int m_infosteps
underlying operator is time stepping
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
int m_nvec
Dimension of Krylov subspace.
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
int m_nits
Number of vectors to test.
Array< OneD, NekDouble > m_imag_evl
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
int m_nfields
interval to dump information if required.
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
static std::string arpackProblemTypeDefault
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DriverArpack(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
static std::string driverLookupId
static std::string arpackProblemTypeLookupIds[]
void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
static std::string className
Name of the class.
static std::string ArpackProblemTypeTrans[]
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
LibUtilities::CommSharedPtr m_comm
Communication object.
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
int m_nequ
number of equations
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.
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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static const NekDouble kNekUnsetDouble
static const NekDouble kNekZeroTol
DriverFactory & GetDriverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::vector< double > z(NPUPPER)
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.