Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | 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 &vComm, const int nscale, const NekSysKey &pKey)
 
 ~NekNonlinSysNewton ()
 
- Public Member Functions inherited from Nektar::LibUtilities::NekNonlinSys
 NekNonlinSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, 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 SetRefResidual (const Array< OneD, const NekDouble > &in)
 
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 ()
 
- Public Member Functions inherited from Nektar::LibUtilities::NekSys
 NekSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, 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)
 
virtual void v_NekSysInitialGuess (const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
 
void SetFlagWarnings (bool in)
 

Static Public Member Functions

static NekNonlinSysSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
 
- Static Public Member Functions inherited from Nektar::LibUtilities::NekNonlinSys
static NekNonlinSysSharedPtr CreateInstance (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, 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 &vComm, const int nDimen, const NekSysKey &pKey)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

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)
 
void CalcInexactNewtonForcing (const int &k, NekDouble &resnormOld, const NekDouble &resnorm, NekDouble &forcing)
 
virtual void v_SetupNekNonlinSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir)
 

Protected Attributes

NekDouble m_SysResNorm0
 
NekDouble m_SysResNorm
 
int 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
 
bool m_ResidualUpdated = false
 
- Protected Attributes inherited from Nektar::LibUtilities::NekSys
int m_maxiter
 Maximum iterations. More...
 
NekDouble m_tolerance
 Tolerance of iterative solver. More...
 
LibUtilities::CommSharedPtr m_Comm
 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...
 

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 vComm,
const int  nscale,
const NekSysKey pKey 
)

Definition at line 54 of file NekNonlinSysNewton.cpp.

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

◆ ~NekNonlinSysNewton()

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

Definition at line 71 of file NekNonlinSysNewton.cpp.

72 {
73 }

Member Function Documentation

◆ CalcInexactNewtonForcing()

void Nektar::LibUtilities::NekNonlinSysNewton::CalcInexactNewtonForcing ( const int &  k,
NekDouble resnormOld,
const NekDouble resnorm,
NekDouble forcing 
)
protected

Definition at line 165 of file NekNonlinSysNewton.cpp.

170 {
171  if (0 == k)
172  {
173  forcing = m_LinSysRelativeTolInNonlin;
174  resnormOld = resnorm;
175  }
176  else
177  {
178  switch(m_InexactNewtonForcing)
179  {
180  case 0:
181  {
182  forcing = m_LinSysRelativeTolInNonlin;
183  break;
184  }
185  case 1:
186  {
187  NekDouble tmpForc = m_forcingGamma *
188  pow((resnorm / resnormOld), m_forcingAlpha);
189  NekDouble tmp = m_forcingGamma *
190  pow(forcing, m_forcingAlpha);
191  if (tmp > 0.1)
192  {
193  forcing = min(m_LinSysRelativeTolInNonlin,
194  max(tmp, tmpForc));
195  }
196  else
197  {
198  forcing = min(m_LinSysRelativeTolInNonlin, tmpForc);
199  }
200 
201  forcing = max(forcing, 1.0E-6);
202  break;
203  }
204  }
205  }
206 }
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 vComm,
const int  nDimen,
const NekSysKey pKey 
)
inlinestatic

Definition at line 54 of file NekNonlinSysNewton.h.

58  {
61  vComm, 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: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 
)
protectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 135 of file NekNonlinSysNewton.cpp.

138 {
139  bool converged = false;
140  NekDouble resratio = 1.0;
141  int ntotal = Residual.size();
142 
143  m_SysResNorm = Vmath::Dot(ntotal, Residual, Residual);
145 
146  if (0 == nIteration)
147  {
149  resratio = 1.0;
150  }
151  else
152  {
153  resratio = m_SysResNorm / m_SysResNorm0;
154  }
155 
157  m_SysResNorm < tol)
158  {
159  converged = true;
160  }
161 
162  return converged;
163 }
LibUtilities::CommSharedPtr m_Comm
Communicate.
Definition: NekSys.h:285
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:1038

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

Referenced by v_SolveSystem().

◆ v_InitObject()

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

Reimplemented from Nektar::LibUtilities::NekNonlinSys.

Definition at line 63 of file NekNonlinSysNewton.cpp.

64 {
66  m_Residual = Array<OneD, NekDouble>(m_SysDimen, 0.0);
67  m_DeltSltn = Array<OneD, NekDouble>(m_SysDimen, 0.0);
68  m_SourceVec = Array<OneD, NekDouble>(m_SysDimen, 0.0);
69 }
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:155
Array< OneD, NekDouble > m_DeltSltn
Definition: NekNonlinSys.h:156
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:157
virtual void v_InitObject()
Definition: NekSys.h:298
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:296

References Nektar::LibUtilities::NekNonlinSys::m_DeltSltn, Nektar::LibUtilities::NekNonlinSys::m_Residual, Nektar::LibUtilities::NekNonlinSys::m_SourceVec, 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 
)
protectedvirtual

Reimplemented from Nektar::LibUtilities::NekNonlinSys.

Definition at line 208 of file NekNonlinSysNewton.cpp.

212 {
213  boost::ignore_unused(nGlobal, nDir);
214 
215  ASSERTL0(0 == nDir, "0 != nDir not tested");
216  ASSERTL0(m_SysDimen == nGlobal, "m_SysDimen!=nGlobal");
217 
218  m_SourceVec = pSource;
219  Vmath::Vcopy(nGlobal-nDir, pSource, 1, m_SourceVec, 1);
220 
221  if (!m_ResidualUpdated)
222  {
224  m_ResidualUpdated = true;
225  }
226  m_linsol->SetSysOperators(m_operator);
227 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:143
NekSysOperators m_operator
Operators.
Definition: NekSys.h:294
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:128
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References ASSERTL0, Nektar::LibUtilities::NekSysOperators::DoNekSysResEval(), Nektar::LibUtilities::NekNonlinSys::m_linsol, Nektar::LibUtilities::NekSys::m_operator, Nektar::LibUtilities::NekNonlinSys::m_Residual, Nektar::LibUtilities::NekNonlinSys::m_ResidualUpdated, Nektar::LibUtilities::NekNonlinSys::m_SourceVec, Nektar::LibUtilities::NekSys::m_SysDimen, and Vmath::Vcopy().

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 
)
protectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 78 of file NekNonlinSysNewton.cpp.

82 {
83  boost::ignore_unused(factor);
84 
85  int nwidthcolm = 11;
86 
87  v_SetupNekNonlinSystem(nGlobal, pInput, pInput, nDir);
88 
89  m_Solution = pOutput;
90  Vmath::Vcopy(nGlobal-nDir, pInput, 1, m_Solution, 1);
91 
92  int ntotal = nGlobal - nDir;
93  m_NtotLinSysIts = 0;
94 
95  int NttlNonlinIte = 0;
96  m_converged = false;
97 
98  NekDouble resnormOld = 0.0;
99  for (int k = 0; k < m_maxiter; ++k)
100  {
102  if (m_converged)
103  break;
104 
105  NekDouble LinSysRelativeIteTol;
106  CalcInexactNewtonForcing(k, resnormOld, m_SysResNorm,
107  LinSysRelativeIteTol);
108 
109  resnormOld = m_SysResNorm;
110 
111  NekDouble LinSysTol = LinSysRelativeIteTol * sqrt(m_SysResNorm);
112  int ntmpLinSysIts =
113  m_linsol->SolveSystem(ntotal, m_Residual, m_DeltSltn, 0, LinSysTol);
114  m_NtotLinSysIts += ntmpLinSysIts;
115  Vmath::Vsub(ntotal, m_Solution, 1, m_DeltSltn, 1, m_Solution, 1);
116  NttlNonlinIte++;
118  }
119 
120  if ( ((!m_converged) || m_verbose) && m_root && m_FlagWarnings)
121  {
123  " # Nonlinear solver not converge in DoImplicitSolve");
124  cout << right << scientific << setw(nwidthcolm)
125  << setprecision(nwidthcolm - 6)
126  << " * Newton-Its converged (RES=" << sqrt(m_SysResNorm)
127  << " Res/(DtRHS): " << sqrt(m_SysResNorm / m_SysResNorm0)
128  << " with " << setw(3) << NttlNonlinIte << " Non-Its)" << endl;
129  }
130 
131  m_ResidualUpdated = false;
132  return NttlNonlinIte;
133 }
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:223
Array< OneD, NekDouble > m_Solution
Definition: NekNonlinSys.h:154
virtual bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol)
void CalcInexactNewtonForcing(const int &k, NekDouble &resnormOld, const NekDouble &resnorm, NekDouble &forcing)
virtual void v_SetupNekNonlinSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir)
bool m_root
Root if parallel.
Definition: NekSys.h:289
bool m_verbose
Verbose.
Definition: NekSys.h:291
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:287
int m_maxiter
Maximum iterations.
Definition: NekSys.h:281
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:372
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

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::NekNonlinSys::m_ResidualUpdated, 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:200
static NekNonlinSysSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, 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 80 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_forcingGamma

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

Definition at line 79 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_InexactNewtonForcing

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

Definition at line 78 of file NekNonlinSysNewton.h.

Referenced by CalcInexactNewtonForcing().

◆ m_SysResNorm

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

Definition at line 76 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().

◆ m_SysResNorm0

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

Definition at line 75 of file NekNonlinSysNewton.h.

Referenced by v_ConvergenceCheck(), and v_SolveSystem().