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 ()
 
- Public Member Functions inherited from Nektar::LibUtilities::NekNonlinSys
 NekNonlinSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
 ~NekNonlinSys ()
 
const Array< OneD, const NekDouble > & GetRefSolution () const
 
const Array< OneD, const NekDouble > & GetRefResidual () const
 
const Array< OneD, const NekDouble > & GetRefSourceVec () const
 
void SetNekNonlinSysTolerance (const NekDouble in)
 
void SetNekNonlinSysMaxIterations (const unsigned int in)
 
const NekLinSysIterSharedPtrGetLinSys ()
 
void SetNonlinIterTolRelativeL2 (const NekDouble in)
 
void SetLinSysRelativeTolInNonlin (const NekDouble in)
 
int GetNtotLinSysIts ()
 
void SetupNekNonlinSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir)
 
- Public Member Functions inherited from Nektar::LibUtilities::NekSys
 NekSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
void InitObject ()
 
virtual ~NekSys ()
 
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

virtual void v_InitObject () override
 
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
 
virtual bool v_ConvergenceCheck (const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol) override
 
virtual void v_SetupNekNonlinSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir) override
 
virtual void v_InitObject ()
 
virtual void v_SetupNekNonlinSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir)=0
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekSys
virtual void v_InitObject ()
 
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))
 
- 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 48 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 54 of file NekNonlinSysNewton.cpp.

58 : NekNonlinSys(pSession, vRowComm, nscale, pKey)
59{
60}
NekNonlinSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)

◆ ~NekNonlinSysNewton()

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

Definition at line 69 of file NekNonlinSysNewton.cpp.

70{
71}

Member Function Documentation

◆ CalcInexactNewtonForcing()

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

Definition at line 162 of file NekNonlinSysNewton.cpp.

164{
165 if (k == 0 || !m_InexactNewtonForcing)
166 {
168 }
169 else
170 {
171 NekDouble tmpForc =
172 m_forcingGamma * pow((resnorm / resnormOld), m_forcingAlpha);
173 return max(min(m_LinSysRelativeTolInNonlin, tmpForc), 1.0E-6);
174 }
175}
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 54 of file NekNonlinSysNewton.h.

58 {
61 pSession, vRowComm, nDimen, pKey);
62 p->InitObject();
63 return p;
64 }
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:48

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 132 of file NekNonlinSysNewton.cpp.

135{
136 bool converged = false;
137 NekDouble resratio = 1.0;
138 int ntotal = Residual.size();
139
140 m_SysResNorm = Vmath::Dot(ntotal, Residual, Residual);
142
143 if (nIteration == 0)
144 {
146 resratio = 1.0;
147 }
148 else
149 {
150 resratio = m_SysResNorm / m_SysResNorm0;
151 }
152
154 m_SysResNorm < tol * tol)
155 {
156 converged = true;
157 }
158
159 return converged;
160}
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:302
T Dot(int n, const T *w, const T *x)
dot (vector times vector): z = w*x
Definition: Vmath.cpp:1095

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 62 of file NekNonlinSysNewton.cpp.

63{
65 m_Residual = Array<OneD, NekDouble>(m_SysDimen, 0.0);
66 m_DeltSltn = Array<OneD, NekDouble>(m_SysDimen, 0.0);
67}
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:133
Array< OneD, NekDouble > m_DeltSltn
Definition: NekNonlinSys.h:134
virtual void v_InitObject()
Definition: NekSys.h:315
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:313

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

◆ v_SetupNekNonlinSystem()

void Nektar::LibUtilities::NekNonlinSysNewton::v_SetupNekNonlinSystem ( const int  nGlobal,
const Array< OneD, const NekDouble > &  pInput,
const Array< OneD, const NekDouble > &  pSource,
const int  nDir 
)
overrideprotectedvirtual

Implements Nektar::LibUtilities::NekNonlinSys.

Definition at line 177 of file NekNonlinSysNewton.cpp.

180{
181 boost::ignore_unused(pInput);
182
183 ASSERTL0(0 == nDir, "0 != nDir not tested");
184 ASSERTL0(m_SysDimen == nGlobal, "m_SysDimen!=nGlobal");
185
186 m_SourceVec = pSource;
187
188 m_linsol->SetSysOperators(m_operator);
189}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:135
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:122
NekSysOperators m_operator
Operators.
Definition: NekSys.h:311

References ASSERTL0, Nektar::LibUtilities::NekNonlinSys::m_linsol, Nektar::LibUtilities::NekSys::m_operator, Nektar::LibUtilities::NekNonlinSys::m_SourceVec, and Nektar::LibUtilities::NekSys::m_SysDimen.

Referenced by v_SolveSystem().

◆ 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 76 of file NekNonlinSysNewton.cpp.

80{
81 boost::ignore_unused(factor);
82
83 v_SetupNekNonlinSystem(nGlobal, pInput, pInput, nDir);
84
85 m_Solution = pOutput;
86 Vmath::Vcopy(nGlobal - nDir, pInput, 1, m_Solution, 1);
87
88 int ntotal = nGlobal - nDir;
90 m_converged = false;
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 NekDouble LinSysTol = LinSysRelativeIteTol * sqrt(m_SysResNorm);
107 resnormOld = m_SysResNorm;
108
109 int ntmpLinSysIts =
110 m_linsol->SolveSystem(ntotal, m_Residual, m_DeltSltn, 0, LinSysTol);
111 m_NtotLinSysIts += ntmpLinSysIts;
112
113 Vmath::Vsub(ntotal, m_Solution, 1, m_DeltSltn, 1, m_Solution, 1);
114 }
115
116 if (((!m_converged) || m_verbose) && m_root && m_FlagWarnings)
117 {
118 int nwidthcolm = 11;
119
121 " # Nonlinear solver not converge in DoImplicitSolve");
122 cout << right << scientific << setw(nwidthcolm)
123 << setprecision(nwidthcolm - 6)
124 << " * Newton-Its converged (RES=" << sqrt(m_SysResNorm)
125 << " Res/(DtRHS): " << sqrt(m_SysResNorm / m_SysResNorm0)
126 << " with " << setw(3) << NttlNonlinIte << " Non-Its)" << endl;
127 }
128
129 return NttlNonlinIte;
130}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:222
Array< OneD, NekDouble > m_Solution
Definition: NekNonlinSys.h:132
virtual void v_SetupNekNonlinSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir) override
NekDouble CalcInexactNewtonForcing(const int &k, const NekDouble &resnormOld, const NekDouble &resnorm)
virtual 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:306
bool m_verbose
Verbose.
Definition: NekSys.h:308
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:304
int m_maxiter
Maximum iterations.
Definition: NekSys.h:298
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:136
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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.cpp:414
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References 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, m_SysResNorm, m_SysResNorm0, Nektar::LibUtilities::NekSys::m_verbose, tinysimd::sqrt(), v_ConvergenceCheck(), v_SetupNekNonlinSystem(), 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:198
static NekNonlinSysSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekNonlinSysFactory & GetNekNonlinSysFactory()

Definition at line 66 of file NekNonlinSysNewton.h.

◆ m_forcingAlpha

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

Definition at line 79 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_forcingGamma

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

Definition at line 78 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_InexactNewtonForcing

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

Definition at line 77 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_SysResNorm

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

Definition at line 75 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().

◆ m_SysResNorm0

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

Definition at line 74 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().