36 #ifndef NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_NekSys_H
37 #define NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_NekSys_H
45 namespace LibUtilities
99 template <
typename FuncPo
interT,
typename ObjectPo
interT>
103 std::bind(
func, obj, std::placeholders::_1, std::placeholders::_2,
104 std::placeholders::_3);
106 template <
typename FuncPo
interT,
typename ObjectPo
interT>
110 std::bind(
func, obj, std::placeholders::_1, std::placeholders::_2,
111 std::placeholders::_3);
113 template <
typename FuncPo
interT,
typename ObjectPo
interT>
117 std::bind(
func, obj, std::placeholders::_1, std::placeholders::_2,
118 std::placeholders::_3);
120 template <
typename FuncPo
interT,
typename ObjectPo
interT>
124 std::bind(
func, obj, std::placeholders::_1, std::placeholders::_2,
125 std::placeholders::_3, std::placeholders::_4);
129 const bool &flag =
false)
const
136 const bool &flag =
false)
const
143 const bool &flag =
false)
const
157 const bool &flag =
false)
const
220 class NekSys :
public std::enable_shared_from_this<NekSys>
260 return v_SolveSystem(nGlobal, pInput, pOutput, nDir, tol, factor);
307 boost::ignore_unused(nGlobal, pInput, pOutput, nDir, tol, factor);
308 ASSERTL0(
false,
"LinSysIterSolver NOT CORRECT.");
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define LIB_UTILITIES_EXPORT
void SetFlagWarnings(bool in)
bool m_root
Root if parallel.
virtual void v_NekSysInitialGuess(const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
void SetSysOperators(const NekSysOperators &in)
virtual void v_InitObject()
bool ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol=1.0E-7)
static NekSysSharedPtr CreateInstance(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
NekDouble m_tolerance
Tolerance of iterative solver.
NekSysOperators m_operator
Operators.
const NekSysOperators & GetSysOperators()
NekSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
int m_SysDimen
The dimension of the system.
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)
virtual bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol)
LibUtilities::CommSharedPtr m_Comm
Communicate.
int SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir, const NekDouble tol=1.0E-7, const NekDouble factor=1.0)
bool m_converged
Whether the iteration has been converged.
int m_maxiter
Maximum iterations.
NekDouble m_NonlinIterTolRelativeL2
int m_NekNonlinSysMaxIterations
NekDouble m_NekLinSysTolerance
bool m_NekLinSysLeftPrecon
NekDouble m_LinSysRelativeTolInNonlin
std::string m_LinSysIterSolverTypeInNonlin
NekDouble m_NekNonlinSysTolerance
bool m_NekLinSysRightPrecon
int m_KrylovMaxHessMatBand
int m_NekLinSysMaxIterations
void DefineNekSysFixPointIte(FuncPointerT func, ObjectPointerT obj)
Array< OneD, FunctorType2 > FunctorType2Array
NekSysOperators & operator=(const NekSysOperators &in)
static const int nfunctor2
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
static const int nfunctor1
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
const Array< OneD, const NekDouble > InArrayType
FunctorType2Array m_functors2
std::function< void(InArrayType &, InArrayType &, OutArrayType &, const bool &)> FunctorType2
Array< OneD, NekDouble > OutArrayType
std::function< void(InArrayType &, OutArrayType &, const bool &)> FunctorType1
NekSysOperators(const NekSysOperators &in)
FunctorType1Array m_functors1
void DoNekSysFixPointIte(InArrayType &rhs, InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Array< OneD, FunctorType1 > FunctorType1Array
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< NekSys > NekSysSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
static const NekDouble kNekIterativeTol
The above copyright notice and this permission notice shall be included.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)