43 namespace LibUtilities
50 string NekLinSysIterGMRES::className =
52 "GMRES", NekLinSysIterGMRES::create,
"NekLinSysIterGMRES solver.");
54 NekLinSysIterGMRES::NekLinSysIterGMRES(
60 std::vector<std::string> variables(1);
61 variables[0] = pSession->GetVariable(0);
62 string variable = variables[0];
69 if (pSession->DefinesGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand"))
72 pSession->GetGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand")
84 int GMRESCentralDifference = 0;
85 pSession->LoadParameter(
"GMRESCentralDifference", GMRESCentralDifference,
88 switch (GMRESCentralDifference)
123 boost::ignore_unused(tol);
127 int niterations =
DoGMRES(nGlobal, pInput, pOutput, nDir);
159 int nNonDir = nGlobal - nDir;
171 bool restarted =
false;
172 bool truncted =
false;
181 for (
int nrestart = 0; nrestart <
m_maxrestart; ++nrestart)
184 DoGmresRestart(restarted, truncted, nGlobal, pInput, pOutput, nDir);
197 Vmath::Svtvp(nNonDir, -1.0, &r0[0] + nDir, 1, &pInput[0] + nDir, 1,
206 cout << std::scientific << std::setw(nwidthcolm)
207 << std::setprecision(nwidthcolm - 8)
213 cout <<
" WITH (GMRES eps = " << eps <<
" REAL eps= " << eps1
218 cout <<
" CONVERGED" << endl;
222 cout <<
" WARNING: Exceeded maxIt" << endl;
235 const bool restarted,
const bool truncted,
const int nGlobal,
239 int nNonDir = nGlobal - nDir;
298 Vmath::Svtvp(nNonDir, beta, &r0[0] + nDir, 1, &pInput[0] + nDir, 1,
304 Vmath::Vcopy(nNonDir, &pInput[0] + nDir, 1, &r0[0] + nDir, 1);
328 vExchange =
Vmath::Dot2(nNonDir, &pInput[0] + nDir,
329 &pInput[0] + nDir, &
m_map[0] + nDir);
351 if (truncted && (starttem) > 0)
353 id_start[nd] = starttem;
362 alpha = 1.0 / eta[0];
364 Vmath::Smul(nNonDir, alpha, &r0[0] + nDir, 1, &V_total[0][0] + nDir, 1);
375 Vsingle2 = V_total[nd + 1];
380 tmp1 = V_total[nd] + nDir;
381 tmp2 = Vsingle1 + nDir;
386 Vsingle1 = V_total[nd];
390 starttem = id_start[idtem];
391 endtem = id_end[idtem];
393 DoArnoldi(starttem, endtem, nGlobal, nDir, V_total, Vsingle1, Vsingle2,
398 starttem = starttem - 1;
401 hsingle2 = Upper[nd];
406 eps = eta[nd + 1] * eta[nd + 1];
436 for (
int i = 0; i < nswp; ++i)
440 &solution[0] + nDir, 1, &solution[0] + nDir, 1);
445 tmp1 = solution + nDir;
446 tmp2 = solution + nDir;
449 Vmath::Vadd(nNonDir, &solution[0] + nDir, 1, &pOutput[0] + nDir, 1,
450 &pOutput[0] + nDir, 1);
457 const int nGlobal,
const int nDir,
472 int nNonDir = nGlobal - nDir;
492 numbertem = starttem;
493 for (
int i = starttem; i < endtem; ++i)
496 Vmath::Dot2(nNonDir, &w[0] + nDir, &V_local[numbertem][0] + nDir,
500 hsingle[i] = vExchange;
502 beta = -1.0 * vExchange;
503 Vmath::Svtvp(nNonDir, beta, &V_local[numbertem][0] + nDir, 1,
504 &w[0] + nDir, 1, &w[0] + nDir, 1);
505 numbertem = numbertem + 1;
515 hsingle[endtem] =
sqrt(vExchange);
517 alpha = 1.0 / hsingle[endtem];
518 Vmath::Smul(nNonDir, alpha, &w[0] + nDir, 1, &Vsingle2[0] + nDir, 1);
523 const int nGlobal,
const int nDir,
529 boost::ignore_unused(nGlobal, nDir);
533 int idtem = endtem - 1;
539 for (
int i = starttem; i < idtem; ++i)
541 temp_dbl = c[i] * hsingle[i] - s[i] * hsingle[i + 1];
542 hsingle[i + 1] = s[i] * hsingle[i] + c[i] * hsingle[i + 1];
543 hsingle[i] = temp_dbl;
546 hh = hsingle[endtem];
552 else if (
abs(hh) >
abs(dd))
555 s[idtem] = 1.0 /
sqrt(1.0 + temp_dbl * temp_dbl);
556 c[idtem] = temp_dbl * s[idtem];
561 c[idtem] = 1.0 /
sqrt(1.0 + temp_dbl * temp_dbl);
562 s[idtem] = temp_dbl * c[idtem];
565 hsingle[idtem] = c[idtem] * hsingle[idtem] - s[idtem] * hsingle[endtem];
566 hsingle[endtem] = 0.0;
568 temp_dbl = c[idtem] * eta[idtem] - s[idtem] * eta[endtem];
569 eta[endtem] = s[idtem] * eta[idtem] + c[idtem] * eta[endtem];
570 eta[idtem] = temp_dbl;
583 int maxid = number - 1;
585 y[maxid] = b[maxid] /
A[maxid][maxid];
587 for (
int i = maxid - 1; i > -1; --i)
591 for (
int j = i + 1; j < number; ++j)
594 sum = sum - y[j] *
A[j][i];
596 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.
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)
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)
int m_KrylovMaxHessMatBand
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)
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()
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
NekDouble m_rhs_magnitude
Dot product of rhs to normalise stopping criterion.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
virtual void v_InitObject()
Array< OneD, int > m_map
Global to universal unique map.
bool m_root
Root if parallel.
NekDouble m_tolerance
Tolerance of iterative solver.
NekSysOperators m_operator
Operators.
LibUtilities::CommSharedPtr m_Comm
Communicate.
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
static const NekDouble kNekUnsetDouble
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)
vvtvp (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)