60 std::vector<std::string> variables(1);
61 variables[0] = pSession->GetVariable(0);
62 string variable = variables[0];
66 pSession->MatchSolverInfo(
"GMRESRightPrecon",
"False",
70 if (pSession->DefinesGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand"))
73 pSession->GetGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand")
85 int GMRESCentralDifference = 0;
86 pSession->LoadParameter(
"GMRESCentralDifference", GMRESCentralDifference,
89 switch (GMRESCentralDifference)
126 int niterations =
DoGMRES(nGlobal, pInput, pOutput, nDir);
153 int nNonDir = nGlobal - nDir;
165 bool restarted =
false;
166 bool truncted =
false;
173 for (
int nrestart = 0; nrestart <
m_maxrestart; ++nrestart)
176 DoGmresRestart(restarted, truncted, nGlobal, pInput, pOutput, nDir);
189 Vmath::Svtvp(nNonDir, -1.0, &r0[0] + nDir, 1, &pInput[0] + nDir, 1,
200 cout << std::scientific << std::setw(nwidthcolm)
201 << std::setprecision(nwidthcolm - 8)
207 cout <<
" WITH (GMRES eps = " << eps <<
" REAL eps= " << eps1
212 cout <<
" CONVERGED" << endl;
216 cout <<
" WARNING: Exceeded maxIt" << endl;
229 const bool restarted,
const bool truncted,
const int nGlobal,
233 int nNonDir = nGlobal - nDir;
287 Vmath::Vcopy(nNonDir, &pInput[0] + nDir, 1, &r0[0] + nDir, 1);
322 vExchange =
Vmath::Dot2(nNonDir, &pInput[0] + nDir,
323 &pInput[0] + nDir, &
m_map[0] + nDir);
344 if (truncted && (starttem) > 0)
346 id_start[nd] = starttem;
355 alpha = 1.0 / eta[0];
359 Vmath::Smul(nNonDir, alpha, &r0[0] + nDir, 1, &V_total[0][0] + nDir, 1);
373 Vsingle2 = V_total[nd + 1];
378 tmp1 = V_total[nd] + nDir;
379 tmp2 = Vsingle1 + nDir;
384 Vsingle1 = V_total[nd];
389 starttem = id_start[idtem];
390 endtem = id_end[idtem];
392 DoArnoldi(starttem, endtem, nGlobal, nDir, V_total, Vsingle1, Vsingle2,
397 starttem = starttem - 1;
400 hsingle2 = Upper[nd];
405 eps = eta[nd + 1] * eta[nd + 1];
435 for (
int i = 0; i < nswp; ++i)
437 Vmath::Svtvp(nNonDir, y_total[i], &V_total[i][0] + nDir, 1,
438 solution.get(), 1, solution.get(), 1);
447 Vmath::Vadd(nNonDir, solution.get(), 1, &pOutput[0] + nDir, 1,
448 &pOutput[0] + nDir, 1);
455 const int nGlobal,
const int nDir,
470 int nNonDir = nGlobal - nDir;
494 numbertem = starttem;
495 for (
int i = starttem; i < endtem; ++i)
498 Vmath::Dot2(nNonDir, &
w[0] + nDir, &V_local[numbertem][0] + nDir,
502 hsingle[i] = vExchange;
504 beta = -1.0 * vExchange;
506 &
w[0] + nDir, 1, &
w[0] + nDir, 1);
507 numbertem = numbertem + 1;
517 hsingle[endtem] =
sqrt(vExchange);
519 alpha = 1.0 / hsingle[endtem];
520 Vmath::Smul(nNonDir, alpha, &
w[0] + nDir, 1, &Vsingle2[0] + nDir, 1);
525 const int nGlobal,
const int nDir,
531 boost::ignore_unused(nGlobal, nDir);
535 int idtem = endtem - 1;
541 for (
int i = starttem; i < idtem; ++i)
543 temp_dbl = c[i] * hsingle[i] - s[i] * hsingle[i + 1];
544 hsingle[i + 1] = s[i] * hsingle[i] + c[i] * hsingle[i + 1];
545 hsingle[i] = temp_dbl;
548 hh = hsingle[endtem];
554 else if (
abs(hh) >
abs(dd))
557 s[idtem] = 1.0 /
sqrt(1.0 + temp_dbl * temp_dbl);
558 c[idtem] = temp_dbl * s[idtem];
563 c[idtem] = 1.0 /
sqrt(1.0 + temp_dbl * temp_dbl);
564 s[idtem] = temp_dbl * c[idtem];
567 hsingle[idtem] = c[idtem] * hsingle[idtem] - s[idtem] * hsingle[endtem];
568 hsingle[endtem] = 0.0;
570 temp_dbl = c[idtem] * eta[idtem] - s[idtem] * eta[endtem];
571 eta[endtem] = s[idtem] * eta[idtem] + c[idtem] * eta[endtem];
572 eta[idtem] = temp_dbl;
585 int maxid = number - 1;
587 y[maxid] = b[maxid] /
A[maxid][maxid];
589 for (
int i = maxid - 1; i > -1; --i)
593 for (
int j = i + 1; j < number; ++j)
596 sum = sum - y[j] *
A[j][i];
598 y[i] = sum /
A[i][i];
#define WARNINGL1(condition, msg)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
NekLinSysIterGMRES(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey=NekSysKey())
int DoGMRES(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
Actual iterative solve-GMRES.
void DoBackward(const int number, Array< OneD, Array< OneD, NekDouble > > &A, const Array< OneD, const NekDouble > &b, Array< OneD, NekDouble > &y)
void DoArnoldi(const int starttem, const int endtem, const int nGlobal, const int nDir, Array< OneD, Array< OneD, NekDouble > > &V_local, Array< OneD, NekDouble > &Vsingle1, Array< OneD, NekDouble > &Vsingle2, Array< OneD, NekDouble > &hsingle)
virtual int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir, const NekDouble tol, const NekDouble factor) override
static std::string className
int m_KrylovMaxHessMatBand
void DoGivensRotation(const int starttem, const int endtem, const int nGlobal, const int nDir, Array< OneD, NekDouble > &c, Array< OneD, NekDouble > &s, Array< OneD, NekDouble > &hsingle, Array< OneD, NekDouble > &eta)
virtual void v_InitObject() override
bool m_NekLinSysRightPrecon
NekDouble DoGmresRestart(const bool restarted, const bool truncted, const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir)
Actual iterative gmres solver for one restart.
bool m_NekLinSysLeftPrecon
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
virtual void v_InitObject() override
NekDouble m_rhs_magnitude
Dot product of rhs to normalise stopping criterion.
Array< OneD, int > m_map
Global to universal unique map.
bool m_root
Root if parallel.
LibUtilities::CommSharedPtr m_rowComm
Communicate.
NekDouble m_tolerance
Tolerance of iterative solver.
NekSysOperators m_operator
Operators.
bool m_converged
Whether the iteration has been converged.
int m_maxiter
Maximum iterations.
bool m_NekLinSysLeftPrecon
bool m_NekLinSysRightPrecon
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
@ beta
Gauss Radau pinned at x=-1,.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
static const NekDouble kNekUnsetDouble
std::vector< double > w(NPUPPER)
The above copyright notice and this permission notice shall be included.
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
T Dot2(int n, const T *w, const T *x, const int *y)
dot2 (vector times vector times vector): z = w*x*y
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > abs(scalarT< T > in)
scalarT< T > sqrt(scalarT< T > in)