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 &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 ()
 
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 &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)
 
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::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 () 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
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekSys
virtual void v_NekSysInitialGuess (const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
 

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...
 

Private Member Functions

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

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 }
NekNonlinSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)

◆ ~NekNonlinSysNewton()

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

Definition at line 70 of file NekNonlinSysNewton.cpp.

71 {
72 }

Member Function Documentation

◆ CalcInexactNewtonForcing()

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

Definition at line 164 of file NekNonlinSysNewton.cpp.

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

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 134 of file NekNonlinSysNewton.cpp.

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

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 ( )
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  m_SourceVec = Array<OneD, NekDouble>(m_SysDimen, 0.0);
68 }
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:143
Array< OneD, NekDouble > m_DeltSltn
Definition: NekNonlinSys.h:144
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:145
virtual void v_InitObject()
Definition: NekSys.h:293
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:291

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

Implements Nektar::LibUtilities::NekNonlinSys.

Definition at line 205 of file NekNonlinSysNewton.cpp.

208 {
209  boost::ignore_unused(nGlobal, nDir);
210 
211  ASSERTL0(0 == nDir, "0 != nDir not tested");
212  ASSERTL0(m_SysDimen == nGlobal, "m_SysDimen!=nGlobal");
213 
214  m_SourceVec = pSource;
215  Vmath::Vcopy(nGlobal - nDir, pSource, 1, m_SourceVec, 1);
216 
217  if (!m_ResidualUpdated)
218  {
220  m_ResidualUpdated = true;
221  }
222  m_linsol->SetSysOperators(m_operator);
223 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:132
NekSysOperators m_operator
Operators.
Definition: NekSys.h:289
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:1255

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

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 77 of file NekNonlinSysNewton.cpp.

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

int 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().