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

#include <NekNonlinSys.h>

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

Public Member Functions

 NekNonlinSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
 
 ~NekNonlinSys ()
 
virtual void v_SetupNekNonlinSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, const Array< OneD, const NekDouble > &pSource, const int nDir)
 
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 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)
 

Protected Member Functions

virtual void v_InitObject ()
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekSys
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)
 

Protected Attributes

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

Detailed Description

Solves a nonlinear system using iterative methods.

Definition at line 57 of file NekNonlinSys.h.

Constructor & Destructor Documentation

◆ NekNonlinSys()

Nektar::LibUtilities::NekNonlinSys::NekNonlinSys ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vComm,
const int  nDimen,
const NekSysKey pKey 
)

Definition at line 55 of file NekNonlinSys.cpp.

59  : NekSys(pSession, vComm, nDimen, pKey)
60 {
61  std::vector<std::string> variables(1);
62  variables[0] = pSession->GetVariable(0);
63  string variable = variables[0];
64 
65  if (pSession->DefinesGlobalSysSolnInfo(variable, "NekNonlinSysTolerance"))
66  {
67  m_tolerance = boost::lexical_cast<NekDouble>(
68  pSession->GetGlobalSysSolnInfo(variable, "NekNonlinSysTolerance")
69  .c_str());
70  }
71  else
72  {
73  pSession->LoadParameter("NekNonlinSysTolerance", m_tolerance,
74  pKey.m_NekNonlinSysTolerance);
75  }
76 
77  if (pSession->DefinesGlobalSysSolnInfo(variable,
78  "NekNonlinSysMaxIterations"))
79  {
80  m_maxiter = boost::lexical_cast<int>(
81  pSession
82  ->GetGlobalSysSolnInfo(variable, "NekNonlinSysMaxIterations")
83  .c_str());
84  }
85  else
86  {
87  pSession->LoadParameter("NekNonlinSysMaxIterations", m_maxiter,
88  pKey.m_NekNonlinSysMaxIterations);
89  }
90 
91  if (pSession->DefinesGlobalSysSolnInfo(variable, "NonlinIterTolRelativeL2"))
92  {
93  m_NonlinIterTolRelativeL2 = boost::lexical_cast<int>(
94  pSession->GetGlobalSysSolnInfo(variable, "NonlinIterTolRelativeL2")
95  .c_str());
96  }
97  else
98  {
99  pSession->LoadParameter("NonlinIterTolRelativeL2",
100  m_NonlinIterTolRelativeL2, pKey.m_NonlinIterTolRelativeL2);
101  }
102 
103  if (pSession->DefinesGlobalSysSolnInfo(variable,
104  "LinSysRelativeTolInNonlin"))
105  {
106  m_LinSysRelativeTolInNonlin = boost::lexical_cast<int>(
107  pSession
108  ->GetGlobalSysSolnInfo(variable, "LinSysRelativeTolInNonlin")
109  .c_str());
110  }
111  else
112  {
113  pSession->LoadParameter("LinSysRelativeTolInNonlin",
115  pKey.m_LinSysRelativeTolInNonlin);
116  }
117 
118  // cout << " m_LinSysRelativeTolInNonlin = " << m_LinSysRelativeTolInNonlin << endl;
119 
120  m_LinSysIterSolverType = pKey.m_LinSysIterSolverTypeInNonlin;
121  if (pSession->DefinesGlobalSysSolnInfo(variable,
122  "LinSysIterSolverTypeInNonlin"))
123  {
125  pSession->GetGlobalSysSolnInfo(variable,
126  "LinSysIterSolverTypeInNonlin");
127  }
128  else
129  {
130  if (pSession->DefinesSolverInfo("LinSysIterSolverTypeInNonlin"))
131  {
133  pSession->GetSolverInfo("LinSysIterSolverTypeInNonlin");
134  }
135  }
136 
139  "NekLinSysIter '" + m_LinSysIterSolverType +
140  "' is not defined.\n");
141 
143  m_LinSysIterSolverType, pSession, m_Comm, m_SysDimen, pKey);
144  m_linsol->SetFlagWarnings(false);
145 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:143
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:283
NekSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm, const int nDimen, const NekSysKey &pKey)
Definition: NekSys.cpp:50
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:296
LibUtilities::CommSharedPtr m_Comm
Communicate.
Definition: NekSys.h:285
int m_maxiter
Maximum iterations.
Definition: NekSys.h:281
NekLinSysIterFactory & GetNekLinSysIterFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::GetNekLinSysIterFactory(), Nektar::LibUtilities::NekSys::m_Comm, m_linsol, m_LinSysIterSolverType, Nektar::LibUtilities::NekSysKey::m_LinSysIterSolverTypeInNonlin, m_LinSysRelativeTolInNonlin, Nektar::LibUtilities::NekSysKey::m_LinSysRelativeTolInNonlin, Nektar::LibUtilities::NekSys::m_maxiter, Nektar::LibUtilities::NekSysKey::m_NekNonlinSysMaxIterations, Nektar::LibUtilities::NekSysKey::m_NekNonlinSysTolerance, m_NonlinIterTolRelativeL2, Nektar::LibUtilities::NekSysKey::m_NonlinIterTolRelativeL2, Nektar::LibUtilities::NekSys::m_SysDimen, and Nektar::LibUtilities::NekSys::m_tolerance.

◆ ~NekNonlinSys()

Nektar::LibUtilities::NekNonlinSys::~NekNonlinSys ( )

Definition at line 152 of file NekNonlinSys.cpp.

153 {
154 }

Member Function Documentation

◆ CreateInstance()

static NekNonlinSysSharedPtr Nektar::LibUtilities::NekNonlinSys::CreateInstance ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vComm,
const int  nDimen,
const NekSysKey pKey 
)
inlinestatic

Definition at line 61 of file NekNonlinSys.h.

65  {
68  nDimen, pKey);
69  return p;
70  }
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.

◆ GetLinSys()

const NekLinSysIterSharedPtr& Nektar::LibUtilities::NekNonlinSys::GetLinSys ( )
inline

Definition at line 122 of file NekNonlinSys.h.

123  {
124  return m_linsol;
125  }

References m_linsol.

◆ GetNtotLinSysIts()

int Nektar::LibUtilities::NekNonlinSys::GetNtotLinSysIts ( )
inline

Definition at line 137 of file NekNonlinSys.h.

138  {
139  return m_NtotLinSysIts;
140  }

References m_NtotLinSysIts.

◆ GetRefResidual()

const Array<OneD, const NekDouble>& Nektar::LibUtilities::NekNonlinSys::GetRefResidual ( ) const
inline

Definition at line 88 of file NekNonlinSys.h.

90  {
91  return m_Residual;
92  }
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:155

References m_Residual.

◆ GetRefSolution()

const Array<OneD, const NekDouble>& Nektar::LibUtilities::NekNonlinSys::GetRefSolution ( ) const
inline

Definition at line 82 of file NekNonlinSys.h.

84  {
85  return m_Solution;
86  }
Array< OneD, NekDouble > m_Solution
Definition: NekNonlinSys.h:154

References m_Solution.

◆ GetRefSourceVec()

const Array<OneD, const NekDouble>& Nektar::LibUtilities::NekNonlinSys::GetRefSourceVec ( ) const
inline

Definition at line 94 of file NekNonlinSys.h.

96  {
97  return m_SourceVec;
98  }
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:157

References m_SourceVec.

◆ SetLinSysRelativeTolInNonlin()

void Nektar::LibUtilities::NekNonlinSys::SetLinSysRelativeTolInNonlin ( const NekDouble  in)
inline

Definition at line 132 of file NekNonlinSys.h.

133  {
135  }

References m_LinSysRelativeTolInNonlin.

◆ SetNekNonlinSysMaxIterations()

void Nektar::LibUtilities::NekNonlinSys::SetNekNonlinSysMaxIterations ( const unsigned int  in)
inline

Definition at line 115 of file NekNonlinSys.h.

117  {
118  m_maxiter = in;
119  }

References Nektar::LibUtilities::NekSys::m_maxiter.

◆ SetNekNonlinSysTolerance()

void Nektar::LibUtilities::NekNonlinSys::SetNekNonlinSysTolerance ( const NekDouble  in)
inline

Definition at line 110 of file NekNonlinSys.h.

111  {
112  m_tolerance = in;
113  }

References Nektar::LibUtilities::NekSys::m_tolerance.

◆ SetNonlinIterTolRelativeL2()

void Nektar::LibUtilities::NekNonlinSys::SetNonlinIterTolRelativeL2 ( const NekDouble  in)
inline

Definition at line 127 of file NekNonlinSys.h.

128  {
130  }

References m_NonlinIterTolRelativeL2.

◆ SetRefResidual()

void Nektar::LibUtilities::NekNonlinSys::SetRefResidual ( const Array< OneD, const NekDouble > &  in)
inline

Definition at line 100 of file NekNonlinSys.h.

102  {
103  ASSERTL0(in.size() == m_SysDimen,
104  "SetRefResidual dimension not correct");
105  Vmath::Vcopy(m_SysDimen, in, 1, m_Residual, 1);
106 
107  m_ResidualUpdated = true;
108  }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References ASSERTL0, m_Residual, m_ResidualUpdated, Nektar::LibUtilities::NekSys::m_SysDimen, and Vmath::Vcopy().

◆ v_InitObject()

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

Reimplemented from Nektar::LibUtilities::NekSys.

Reimplemented in Nektar::LibUtilities::NekNonlinSysNewton.

Definition at line 147 of file NekNonlinSys.cpp.

148 {
150 }
virtual void v_InitObject()
Definition: NekSys.h:298

References Nektar::LibUtilities::NekSys::v_InitObject().

◆ v_SetupNekNonlinSystem()

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

Reimplemented in Nektar::LibUtilities::NekNonlinSysNewton.

Definition at line 156 of file NekNonlinSys.cpp.

160 {
161  boost::ignore_unused(nGlobal, pInput, pSource, nDir);
162  NEKERROR(ErrorUtil::efatal, "v_SetupNekNonlinSystem not defined");
163 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References Nektar::ErrorUtil::efatal, and NEKERROR.

Friends And Related Function Documentation

◆ MemoryManager< NekNonlinSys >

friend class MemoryManager< NekNonlinSys >
friend

Definition at line 55 of file NekNonlinSys.h.

Member Data Documentation

◆ m_DeltSltn

Array<OneD, NekDouble> Nektar::LibUtilities::NekNonlinSys::m_DeltSltn
protected

◆ m_linsol

NekLinSysIterSharedPtr Nektar::LibUtilities::NekNonlinSys::m_linsol
protected

◆ m_LinSysIterSolverType

std::string Nektar::LibUtilities::NekNonlinSys::m_LinSysIterSolverType
protected

Definition at line 148 of file NekNonlinSys.h.

Referenced by NekNonlinSys().

◆ m_LinSysRelativeTolInNonlin

NekDouble Nektar::LibUtilities::NekNonlinSys::m_LinSysRelativeTolInNonlin
protected

◆ m_NonlinIterTolRelativeL2

NekDouble Nektar::LibUtilities::NekNonlinSys::m_NonlinIterTolRelativeL2
protected

◆ m_NtotLinSysIts

int Nektar::LibUtilities::NekNonlinSys::m_NtotLinSysIts = 0
protected

◆ m_Residual

Array<OneD, NekDouble> Nektar::LibUtilities::NekNonlinSys::m_Residual
protected

◆ m_ResidualUpdated

bool Nektar::LibUtilities::NekNonlinSys::m_ResidualUpdated = false
protected

◆ m_Solution

Array<OneD, NekDouble> Nektar::LibUtilities::NekNonlinSys::m_Solution
protected

◆ m_SourceVec

Array<OneD, NekDouble> Nektar::LibUtilities::NekNonlinSys::m_SourceVec
protected

◆ m_totalIterations

int Nektar::LibUtilities::NekNonlinSys::m_totalIterations = 0
protected

Definition at line 150 of file NekNonlinSys.h.