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

#include <NekSys.h>

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

Public Member Functions

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

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)
 
virtual void v_NekSysInitialGuess (const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
 

Protected Attributes

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< NekSys >
 Support creation through MemoryManager. More...
 

Detailed Description

Solves a nonlinear or linear system.

Definition at line 219 of file NekSys.h.

Constructor & Destructor Documentation

◆ NekSys()

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

Definition at line 50 of file NekSys.cpp.

53 {
54  m_tolerance = pKey.m_Tolerance;
55  m_verbose = false;
56  m_root = false;
57  m_Comm = vComm;
58 
59  m_FlagWarnings = true;
60 
61  if (0 == m_Comm->GetRank())
62  {
63  m_root = true;
64  }
65  m_verbose = pSession->DefinesCmdLineArgument("verbose");
66 
67  m_converged = false;
68 
69  m_SysDimen = nDimen;
70 }
bool m_root
Root if parallel.
Definition: NekSys.h:284
bool m_verbose
Verbose.
Definition: NekSys.h:286
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:278
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:291
LibUtilities::CommSharedPtr m_Comm
Communicate.
Definition: NekSys.h:280
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:282

References Nektar::LibUtilities::NekSysKey::m_Tolerance.

◆ ~NekSys()

Nektar::LibUtilities::NekSys::~NekSys ( )
virtual

Definition at line 72 of file NekSys.cpp.

73 {
74 }

Member Function Documentation

◆ ConvergenceCheck()

bool Nektar::LibUtilities::NekSys::ConvergenceCheck ( const int  nIteration,
const Array< OneD, const NekDouble > &  Residual,
const NekDouble  tol = 1.0E-7 
)
inline

Definition at line 262 of file NekSys.h.

265  {
266  return v_ConvergenceCheck(nIteration, Residual, tol);
267  }
virtual bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol)
Definition: NekSys.cpp:76

References v_ConvergenceCheck().

Referenced by Nektar::LibUtilities::NekLinSysIterFixedpointJacobi::v_SolveSystem().

◆ CreateInstance()

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

Definition at line 225 of file NekSys.h.

229  {
231  pSession, vComm, nDimen, pKey);
232  return p;
233  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< NekSys > NekSysSharedPtr
Definition: NekSys.h:215

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

◆ GetSysOperators()

const NekSysOperators& Nektar::LibUtilities::NekSys::GetSysOperators ( )
inline

Definition at line 249 of file NekSys.h.

250  {
251  return m_operator;
252  }
NekSysOperators m_operator
Operators.
Definition: NekSys.h:289

References m_operator.

◆ InitObject()

void Nektar::LibUtilities::NekSys::InitObject ( )
inline

Definition at line 238 of file NekSys.h.

239  {
240  v_InitObject();
241  }
virtual void v_InitObject()
Definition: NekSys.h:293

References v_InitObject().

◆ SetFlagWarnings()

void Nektar::LibUtilities::NekSys::SetFlagWarnings ( bool  in)
inline

Definition at line 269 of file NekSys.h.

270  {
271  m_FlagWarnings = in;
272  }

References m_FlagWarnings.

◆ SetSysOperators()

void Nektar::LibUtilities::NekSys::SetSysOperators ( const NekSysOperators in)
inline

Definition at line 244 of file NekSys.h.

245  {
246  m_operator = in;
247  }

References m_operator.

◆ SolveSystem()

int Nektar::LibUtilities::NekSys::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 
)
inline

Definition at line 254 of file NekSys.h.

258  {
259  return v_SolveSystem(nGlobal, pInput, pOutput, nDir, tol, factor);
260  }
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)
Definition: NekSys.h:297

References v_SolveSystem().

◆ v_ConvergenceCheck()

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

Reimplemented in Nektar::LibUtilities::NekNonlinSysNewton.

Definition at line 76 of file NekSys.cpp.

79 {
80  bool converged = false;
81  int ntotal = Residual.size();
82  boost::ignore_unused(nIteration);
83 
84  NekDouble SysResNorm = Vmath::Dot(ntotal, Residual, Residual);
85  m_Comm->AllReduce(SysResNorm, Nektar::LibUtilities::ReduceSum);
86 
87  if (SysResNorm < tol)
88  {
89  converged = true;
90  }
91  return converged;
92 }
double NekDouble
T Dot(int n, const T *w, const T *x)
dot (vector times vector): z = w*x
Definition: Vmath.cpp:1100

References Vmath::Dot(), and Nektar::LibUtilities::ReduceSum.

Referenced by ConvergenceCheck().

◆ v_InitObject()

virtual void Nektar::LibUtilities::NekSys::v_InitObject ( )
inlineprotectedvirtual

◆ v_NekSysInitialGuess()

void Nektar::LibUtilities::NekSys::v_NekSysInitialGuess ( const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pguess 
)
protectedvirtual

Natural guess

Definition at line 97 of file NekSys.cpp.

99 {
100  size_t ndim = pInput.size();
101  if (pguess.size() != ndim)
102  {
103  pguess = Array<OneD, NekDouble>{ndim};
104  }
105 
106  Vmath::Vcopy(ndim, pInput, 1, pguess, 1);
107 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Vmath::Vcopy().

◆ v_SolveSystem()

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

Reimplemented in Nektar::LibUtilities::NekNonlinSysNewton, Nektar::LibUtilities::NekLinSysIterGMRES, Nektar::LibUtilities::NekLinSysIterFixedpointJacobi, and Nektar::LibUtilities::NekLinSysIterCG.

Definition at line 297 of file NekSys.h.

301  {
302  boost::ignore_unused(nGlobal, pInput, pOutput, nDir, tol, factor);
303  ASSERTL0(false, "LinSysIterSolver NOT CORRECT.");
304  return 0;
305  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

References ASSERTL0.

Referenced by SolveSystem().

Friends And Related Function Documentation

◆ MemoryManager< NekSys >

friend class MemoryManager< NekSys >
friend

Support creation through MemoryManager.

Definition at line 212 of file NekSys.h.

Member Data Documentation

◆ m_Comm

LibUtilities::CommSharedPtr Nektar::LibUtilities::NekSys::m_Comm
protected

◆ m_converged

bool Nektar::LibUtilities::NekSys::m_converged
protected

◆ m_FlagWarnings

bool Nektar::LibUtilities::NekSys::m_FlagWarnings
protected

◆ m_maxiter

int Nektar::LibUtilities::NekSys::m_maxiter
protected

◆ m_operator

NekSysOperators Nektar::LibUtilities::NekSys::m_operator
protected

◆ m_root

bool Nektar::LibUtilities::NekSys::m_root
protected

◆ m_SysDimen

int Nektar::LibUtilities::NekSys::m_SysDimen
protected

◆ m_tolerance

NekDouble Nektar::LibUtilities::NekSys::m_tolerance
protected

◆ m_verbose

bool Nektar::LibUtilities::NekSys::m_verbose
protected