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

#include <NekNonlinSysNewton.h>

Inheritance diagram for Nektar::LibUtilities::NekNonlinSysNewton:
[legend]

Public Member Functions

 NekNonlinSysNewton (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nscale, const NekSysKey &pKey)
 
 ~NekNonlinSysNewton () override=default
 
- Public Member Functions inherited from Nektar::LibUtilities::NekNonlinSys
 NekNonlinSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
 ~NekNonlinSys () override=default
 
const Array< OneD, const NekDouble > & GetRefSolution () const
 
const Array< OneD, const NekDouble > & GetRefResidual () const
 
const Array< OneD, const NekDouble > & GetRefSourceVec () const
 
const NekLinSysIterSharedPtrGetLinSys ()
 
int GetNtotLinSysIts ()
 
void SetNekNonlinSysTolerance (const NekDouble in)
 
void SetNekNonlinSysMaxIterations (const unsigned int in)
 
void SetNonlinIterTolRelativeL2 (const NekDouble in)
 
void SetLinSysRelativeTolInNonlin (const NekDouble in)
 
- Public Member Functions inherited from Nektar::LibUtilities::NekSys
 NekSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
virtual ~NekSys ()=default
 
void InitObject ()
 
void SetSysOperators (const NekSysOperators &in)
 
const NekSysOperatorsGetSysOperators ()
 
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 ConvergenceCheck (const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol=1.0E-7)
 
void SetFlagWarnings (bool in)
 

Static Public Member Functions

static NekNonlinSysSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
- Static Public Member Functions inherited from Nektar::LibUtilities::NekSys
static NekSysSharedPtr CreateInstance (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

void v_InitObject () override
 
void v_SetSysOperators (const NekSysOperators &in) override
 
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
 
bool v_ConvergenceCheck (const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol) override
 
void v_InitObject () override
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekSys
virtual void v_InitObject ()
 
virtual void v_SetSysOperators (const NekSysOperators &in)
 
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)
 
virtual void v_NekSysInitialGuess (const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
 

Protected Attributes

NekDouble m_SysResNorm0
 
NekDouble m_SysResNorm
 
bool m_InexactNewtonForcing = 0
 
NekDouble m_forcingGamma = 1.0
 
NekDouble m_forcingAlpha = 0.5 * (1.0 + sqrt(5.0))
 
- Protected Attributes inherited from Nektar::LibUtilities::NekNonlinSys
NekLinSysIterSharedPtr m_linsol
 
NekDouble m_NonlinIterTolRelativeL2
 
NekDouble m_LinSysRelativeTolInNonlin
 
std::string m_LinSysIterSolverType
 
int m_totalIterations = 0
 
int m_NtotLinSysIts = 0
 
Array< OneD, NekDoublem_Solution
 
Array< OneD, NekDoublem_Residual
 
Array< OneD, NekDoublem_DeltSltn
 
Array< OneD, NekDoublem_SourceVec
 
- Protected Attributes inherited from Nektar::LibUtilities::NekSys
int m_maxiter
 Maximum iterations. More...
 
NekDouble m_tolerance
 Tolerance of iterative solver. More...
 
LibUtilities::CommSharedPtr m_rowComm
 Communicate. More...
 
bool m_converged
 Whether the iteration has been converged. More...
 
bool m_root
 Root if parallel. More...
 
bool m_verbose
 Verbose. More...
 
bool m_FlagWarnings
 
NekSysOperators m_operator
 Operators. More...
 
int m_SysDimen
 The dimension of the system. More...
 

Private Member Functions

NekDouble CalcInexactNewtonForcing (const int &k, const NekDouble &resnormOld, const NekDouble &resnorm)
 

Friends

class MemoryManager< NekNonlinSysNewton >
 Constructor for full direct matrix solve. More...
 

Detailed Description

Solves a nonlinear system using iterative methods.

Definition at line 46 of file NekNonlinSysNewton.h.

Constructor & Destructor Documentation

◆ NekNonlinSysNewton()

Nektar::LibUtilities::NekNonlinSysNewton::NekNonlinSysNewton ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vRowComm,
const int  nscale,
const NekSysKey pKey 
)

Definition at line 51 of file NekNonlinSysNewton.cpp.

55 : NekNonlinSys(pSession, vRowComm, nscale, pKey)
56{
57}
NekNonlinSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)

◆ ~NekNonlinSysNewton()

Nektar::LibUtilities::NekNonlinSysNewton::~NekNonlinSysNewton ( )
overridedefault

Member Function Documentation

◆ CalcInexactNewtonForcing()

NekDouble Nektar::LibUtilities::NekNonlinSysNewton::CalcInexactNewtonForcing ( const int &  k,
const NekDouble resnormOld,
const NekDouble resnorm 
)
private

Definition at line 149 of file NekNonlinSysNewton.cpp.

151{
152 if (k == 0 || !m_InexactNewtonForcing)
153 {
155 }
156 else
157 {
158 NekDouble tmpForc =
159 m_forcingGamma * pow((resnorm / resnormOld), m_forcingAlpha);
160 return max(min(m_LinSysRelativeTolInNonlin, tmpForc), 1.0E-6);
161 }
162}
double NekDouble

References m_forcingAlpha, m_forcingGamma, m_InexactNewtonForcing, and Nektar::LibUtilities::NekNonlinSys::m_LinSysRelativeTolInNonlin.

Referenced by v_SolveSystem().

◆ create()

static NekNonlinSysSharedPtr Nektar::LibUtilities::NekNonlinSysNewton::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vRowComm,
const int  nDimen,
const NekSysKey pKey 
)
inlinestatic

Definition at line 52 of file NekNonlinSysNewton.h.

56 {
59 pSession, vRowComm, nDimen, pKey);
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< NekNonlinSys > NekNonlinSysSharedPtr
Definition: NekNonlinSys.h:46

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

◆ v_ConvergenceCheck()

bool Nektar::LibUtilities::NekNonlinSysNewton::v_ConvergenceCheck ( const int  nIteration,
const Array< OneD, const NekDouble > &  Residual,
const NekDouble  tol 
)
overrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 131 of file NekNonlinSysNewton.cpp.

134{
135 m_SysResNorm = Vmath::Dot(Residual.size(), Residual, Residual);
137
138 if (nIteration == 0)
139 {
141 }
142
145
146 return resratio < restol * restol || m_SysResNorm < tol * tol;
147}
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:293
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761

References Vmath::Dot(), Nektar::LibUtilities::NekNonlinSys::m_NonlinIterTolRelativeL2, Nektar::LibUtilities::NekSys::m_rowComm, m_SysResNorm, m_SysResNorm0, and Nektar::LibUtilities::ReduceSum.

Referenced by v_SolveSystem().

◆ v_InitObject()

void Nektar::LibUtilities::NekNonlinSysNewton::v_InitObject ( )
overrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::NekNonlinSys.

Definition at line 59 of file NekNonlinSysNewton.cpp.

60{
62 m_Residual = Array<OneD, NekDouble>(m_SysDimen, 0.0);
63 m_DeltSltn = Array<OneD, NekDouble>(m_SysDimen, 0.0);
64}
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:124
Array< OneD, NekDouble > m_DeltSltn
Definition: NekNonlinSys.h:125
virtual void v_InitObject()
Definition: NekSys.h:306
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:304

References Nektar::LibUtilities::NekNonlinSys::m_DeltSltn, Nektar::LibUtilities::NekNonlinSys::m_Residual, Nektar::LibUtilities::NekSys::m_SysDimen, and Nektar::LibUtilities::NekSys::v_InitObject().

◆ v_SetSysOperators()

void Nektar::LibUtilities::NekNonlinSysNewton::v_SetSysOperators ( const NekSysOperators in)
overrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 66 of file NekNonlinSysNewton.cpp.

67{
69 m_linsol->SetSysOperators(in);
70}
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:113
virtual void v_SetSysOperators(const NekSysOperators &in)
Definition: NekSys.h:310

References Nektar::LibUtilities::NekNonlinSys::m_linsol, and Nektar::LibUtilities::NekSys::v_SetSysOperators().

◆ v_SolveSystem()

int Nektar::LibUtilities::NekNonlinSysNewton::v_SolveSystem ( const int  nGlobal,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const int  nDir,
const NekDouble  tol,
const NekDouble  factor 
)
overrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 75 of file NekNonlinSysNewton.cpp.

79{
80 ASSERTL0(nDir == 0, "nDir != 0 not tested");
81 ASSERTL0(nGlobal == m_SysDimen, "nGlobal != m_SysDimen");
82
83 int ntotal = nGlobal - nDir;
84
85 m_SourceVec = pInput;
86 m_Solution = pOutput;
88 m_converged = false;
89
90 Vmath::Vcopy(ntotal, pInput, 1, m_Solution, 1);
91
92 NekDouble resnormOld = 0.0;
93 int NttlNonlinIte = 0;
94 for (; NttlNonlinIte < m_maxiter; ++NttlNonlinIte)
95 {
97
98 m_converged = v_ConvergenceCheck(NttlNonlinIte, m_Residual, tol);
99 if (m_converged)
100 {
101 break;
102 }
103
104 NekDouble LinSysRelativeIteTol =
105 CalcInexactNewtonForcing(NttlNonlinIte, resnormOld, m_SysResNorm);
106 resnormOld = m_SysResNorm;
107 m_linsol->setRhsMagnitude(m_SysResNorm);
108 int ntmpLinSysIts = m_linsol->SolveSystem(
109 ntotal, m_Residual, m_DeltSltn, 0, LinSysRelativeIteTol);
110 m_NtotLinSysIts += ntmpLinSysIts;
111
112 Vmath::Vsub(ntotal, m_Solution, 1, m_DeltSltn, 1, m_Solution, 1);
113 }
114
115 if (((!m_converged) || m_verbose) && m_root && m_FlagWarnings)
116 {
117 int nwidthcolm = 11;
118
120 " # Nonlinear solver not converge in DoImplicitSolve");
121 cout << right << scientific << setw(nwidthcolm)
122 << setprecision(nwidthcolm - 6)
123 << " * Newton-Its converged (RES=" << sqrt(m_SysResNorm)
124 << " Res/(DtRHS): " << sqrt(m_SysResNorm / m_SysResNorm0)
125 << " with " << setw(3) << NttlNonlinIte << " Non-Its)" << endl;
126 }
127
128 return NttlNonlinIte;
129}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:126
Array< OneD, NekDouble > m_Solution
Definition: NekNonlinSys.h:123
NekDouble CalcInexactNewtonForcing(const int &k, const NekDouble &resnormOld, const NekDouble &resnorm)
bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol) override
bool m_root
Root if parallel.
Definition: NekSys.h:297
bool m_verbose
Verbose.
Definition: NekSys.h:299
NekSysOperators m_operator
Operators.
Definition: NekSys.h:302
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:295
int m_maxiter
Maximum iterations.
Definition: NekSys.h:289
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:134
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, CalcInexactNewtonForcing(), Nektar::LibUtilities::NekSysOperators::DoNekSysResEval(), Nektar::LibUtilities::NekSys::m_converged, Nektar::LibUtilities::NekNonlinSys::m_DeltSltn, Nektar::LibUtilities::NekSys::m_FlagWarnings, Nektar::LibUtilities::NekNonlinSys::m_linsol, Nektar::LibUtilities::NekSys::m_maxiter, Nektar::LibUtilities::NekNonlinSys::m_NtotLinSysIts, Nektar::LibUtilities::NekSys::m_operator, Nektar::LibUtilities::NekNonlinSys::m_Residual, Nektar::LibUtilities::NekSys::m_root, Nektar::LibUtilities::NekNonlinSys::m_Solution, Nektar::LibUtilities::NekNonlinSys::m_SourceVec, Nektar::LibUtilities::NekSys::m_SysDimen, m_SysResNorm, m_SysResNorm0, Nektar::LibUtilities::NekSys::m_verbose, tinysimd::sqrt(), v_ConvergenceCheck(), Vmath::Vcopy(), Vmath::Vsub(), and WARNINGL0.

Friends And Related Function Documentation

◆ MemoryManager< NekNonlinSysNewton >

friend class MemoryManager< NekNonlinSysNewton >
friend

Constructor for full direct matrix solve.

Definition at line 1 of file NekNonlinSysNewton.h.

Member Data Documentation

◆ className

string Nektar::LibUtilities::NekNonlinSysNewton::className
static
Initial value:
=
"Newton", NekNonlinSysNewton::create, "NekNonlinSysNewton solver.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static NekNonlinSysSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekNonlinSysFactory & GetNekNonlinSysFactory()

Definition at line 64 of file NekNonlinSysNewton.h.

◆ m_forcingAlpha

NekDouble Nektar::LibUtilities::NekNonlinSysNewton::m_forcingAlpha = 0.5 * (1.0 + sqrt(5.0))
protected

Definition at line 78 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_forcingGamma

NekDouble Nektar::LibUtilities::NekNonlinSysNewton::m_forcingGamma = 1.0
protected

Definition at line 77 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_InexactNewtonForcing

bool Nektar::LibUtilities::NekNonlinSysNewton::m_InexactNewtonForcing = 0
protected

Definition at line 76 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_SysResNorm

NekDouble Nektar::LibUtilities::NekNonlinSysNewton::m_SysResNorm
protected

Definition at line 74 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().

◆ m_SysResNorm0

NekDouble Nektar::LibUtilities::NekNonlinSysNewton::m_SysResNorm0
protected

Definition at line 73 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().