Nektar++
NekSys.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NekSys.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: NekSys header
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_NekSys_H
37#define NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_NEK_NekSys_H
38
42#include <iomanip>
44{
45
46// =====================================================================
47// ==== DEFINITION OF THE CLASS NekSysOperators
48// ==== Defines operators needed by iterative solver of nonlinear or
49// ==== linear system
50// =====================================================================
52{
53public:
56
57 typedef std::function<void(InArrayType &, OutArrayType &, const bool &)>
59 typedef std::function<void(InArrayType &, InArrayType &, OutArrayType &,
60 const bool &)>
64 static const int nfunctor1 = 4;
65 static const int nfunctor2 = 1;
66
68 {
69 }
70
73 {
74 for (int i = 0; i < nfunctor1; ++i)
75 {
76 m_functors1[i] = in.m_functors1[i];
77 }
78 for (int i = 0; i < nfunctor2; ++i)
79 {
80 m_functors2[i] = in.m_functors2[i];
81 }
82 }
83
85 {
86 for (int i = 0; i < nfunctor1; ++i)
87 {
88 m_functors1[i] = in.m_functors1[i];
89 }
90 for (int i = 0; i < nfunctor2; ++i)
91 {
92 m_functors2[i] = in.m_functors2[i];
93 }
94
95 return *this;
96 }
97
98 template <typename FuncPointerT, typename ObjectPointerT>
99 void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
100 {
101 m_functors1[0] =
102 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
103 std::placeholders::_3);
104 }
105 template <typename FuncPointerT, typename ObjectPointerT>
106 void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
107 {
108 m_functors1[1] =
109 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
110 std::placeholders::_3);
111 }
112 template <typename FuncPointerT, typename ObjectPointerT>
113 void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
114 {
115 m_functors1[2] =
116 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
117 std::placeholders::_3);
118 }
119 template <typename FuncPointerT, typename ObjectPointerT>
120 void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
121 {
122 m_functors1[3] =
123 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
124 std::placeholders::_3);
125 }
126 template <typename FuncPointerT, typename ObjectPointerT>
127 void DefineNekSysFixPointIte(FuncPointerT func, ObjectPointerT obj)
128 {
129 m_functors2[0] =
130 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
131 std::placeholders::_3, std::placeholders::_4);
132 }
133
134 inline void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray,
135 const bool &flag = false) const
136 {
137 ASSERTL1(m_functors1[0], "DoNekSysResEval should be defined");
138 m_functors1[0](inarray, outarray, flag);
139 }
140
141 inline void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray,
142 const bool &flag = false) const
143 {
144 ASSERTL1(m_functors1[1], "DoNekSysLhsEval should be defined");
145 m_functors1[1](inarray, outarray, flag);
146 }
147
148 inline void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray,
149 const bool &flag = false) const
150 {
151 if (m_functors1[2])
152 {
153 m_functors1[2](inarray, outarray, flag);
154 }
155 else
156 {
157 Vmath::Vcopy(outarray.size(), inarray, 1, outarray, 1);
158 }
159 }
160
162 const bool &flag = false) const
163 {
164 ASSERTL1(m_functors1[3], "DoAssembleLoc should be defined");
165 m_functors1[3](xn, xn1, flag);
166 }
167
169 OutArrayType &xn1,
170 const bool &flag = false) const
171 {
172 ASSERTL1(m_functors2[0], "DoNekSysFixPointIte should be defined");
173 m_functors2[0](rhs, xn, xn1, flag);
174 }
175
176protected:
177 /* Defines operators
178 DoNekSysResEval :
179 evaluations the residual of the Nonlinear/Linear system
180 ie. the residual b-Ax and N(x) for linear and
181 nonlinear systems, respectively
182 May not be used for linear system.
183 DoNekSysLhsEval :
184 evaluations the LHS of the Nonlinear/Linear system (Ax),
185 where A is the matrix and x is solution vector.
186 For linear system A is the coefficient matrix;
187 For nonlinear system A is the coefficient matrix in
188 each nonlinear iterations, for example A is the
189 Jacobian matrix for Newton method;
190 DoNekSysPrecon :
191 Preconditioning operator of the system.
192 DoAssembleLoc :
193 Operator to assemble array and scater it back to the lcoal storage
194 if flag is true will zero Dirichlet conditions
195 DoNekSysFixPointIte :
196 */
198
199 /* Defines one operator
200 DoNekSysFixPointIte :
201 Operator to calculate RHS of fixed point iterations
202 (x^{n+1}=M^{-1}(b-N*x^{n}), with M+N=A).
203 */
205};
206
208{
209public:
210 NekSysKey() = default;
211 ~NekSysKey() = default;
212
223 bool m_DifferenceFlag0 = false;
224 bool m_DifferenceFlag1 = false;
225 std::string m_LinSysIterSolverTypeInNonlin = "GMRES";
226};
227
228class NekSys;
229
230typedef std::shared_ptr<NekSys> NekSysSharedPtr;
231
232class NekSys : public std::enable_shared_from_this<NekSys>
233{
234public:
235 /// Support creation through MemoryManager.
236 friend class MemoryManager<NekSys>;
237
240 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
241 const NekSysKey &pKey)
242 {
244 pSession, vRowComm, nDimen, pKey);
245 return p;
246 }
249 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
250 const NekSysKey &pKey);
251 LIB_UTILITIES_EXPORT virtual ~NekSys() = default;
253 {
254 v_InitObject();
255 }
256
258 {
260 }
261
263 {
264 return m_operator;
265 }
266
268 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
269 Array<OneD, NekDouble> &pOutput, const int nDir,
270 const NekDouble tol = 1.0E-7, const NekDouble factor = 1.0)
271 {
272 return v_SolveSystem(nGlobal, pInput, pOutput, nDir, tol, factor);
273 }
274
276 const int nIteration, const Array<OneD, const NekDouble> &Residual,
277 const NekDouble tol = 1.0E-7)
278 {
279 return v_ConvergenceCheck(nIteration, Residual, tol);
280 }
281
283 {
284 m_FlagWarnings = in;
285 }
286
287protected:
288 /// Maximum iterations
290 /// Tolerance of iterative solver.
292 /// Communicate.
294 /// Whether the iteration has been converged
296 /// Root if parallel
297 bool m_root;
298 /// Verbose
301 /// Operators
303 /// The dimension of the system
305
306 virtual void v_InitObject()
307 {
308 }
309
310 virtual void v_SetSysOperators(const NekSysOperators &in)
311 {
312 m_operator = in;
313 }
314
315 virtual int v_SolveSystem(
316 [[maybe_unused]] const int nGlobal,
317 [[maybe_unused]] const Array<OneD, const NekDouble> &pInput,
318 [[maybe_unused]] Array<OneD, NekDouble> &pOutput,
319 [[maybe_unused]] const int nDir, [[maybe_unused]] const NekDouble tol,
320 [[maybe_unused]] const NekDouble factor)
321 {
322 ASSERTL0(false, "LinSysIterSolver NOT CORRECT.");
323 return 0;
324 }
325
326 virtual bool v_ConvergenceCheck(
327 const int nIteration, const Array<OneD, const NekDouble> &Residual,
328 const NekDouble tol);
329
331 const Array<OneD, const NekDouble> &pInput,
332 Array<OneD, NekDouble> &pguess);
333};
334} // namespace Nektar::LibUtilities
335#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define LIB_UTILITIES_EXPORT
void SetFlagWarnings(bool in)
Definition: NekSys.h:282
bool m_root
Root if parallel.
Definition: NekSys.h:297
const NekSysOperators & GetSysOperators()
Definition: NekSys.h:262
static NekSysSharedPtr CreateInstance(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Definition: NekSys.h:238
virtual void v_NekSysInitialGuess(const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
Definition: NekSys.cpp:73
void SetSysOperators(const NekSysOperators &in)
Definition: NekSys.h:257
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:293
virtual void v_InitObject()
Definition: NekSys.h:306
bool ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol=1.0E-7)
Definition: NekSys.h:275
bool m_verbose
Verbose.
Definition: NekSys.h:299
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:291
NekSysOperators m_operator
Operators.
Definition: NekSys.h:302
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:304
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:315
virtual ~NekSys()=default
NekSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Definition: NekSys.cpp:48
virtual bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol)
Definition: NekSys.cpp:60
virtual void v_SetSysOperators(const NekSysOperators &in)
Definition: NekSys.h:310
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)
Definition: NekSys.h:267
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:295
int m_maxiter
Maximum iterations.
Definition: NekSys.h:289
NekDouble m_NonlinIterTolRelativeL2
Definition: NekSys.h:217
NekDouble m_LinSysRelativeTolInNonlin
Definition: NekSys.h:218
std::string m_LinSysIterSolverTypeInNonlin
Definition: NekSys.h:225
void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:120
void DefineNekSysFixPointIte(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:127
Array< OneD, FunctorType2 > FunctorType2Array
Definition: NekSys.h:63
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:134
void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:99
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:106
const Array< OneD, const NekDouble > InArrayType
Definition: NekSys.h:54
NekSysOperators & operator=(const NekSysOperators &in)
Definition: NekSys.h:84
std::function< void(InArrayType &, InArrayType &, OutArrayType &, const bool &)> FunctorType2
Definition: NekSys.h:61
Array< OneD, NekDouble > OutArrayType
Definition: NekSys.h:55
std::function< void(InArrayType &, OutArrayType &, const bool &)> FunctorType1
Definition: NekSys.h:58
NekSysOperators(const NekSysOperators &in)
Definition: NekSys.h:71
void DoNekSysFixPointIte(InArrayType &rhs, InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:168
void DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:161
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:148
Array< OneD, FunctorType1 > FunctorType1Array
Definition: NekSys.h:62
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:113
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:141
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< NekSys > NekSysSharedPtr
Definition: NekSys.h:230
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
static const NekDouble kNekIterativeTol
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825