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>
43namespace Nektar
44{
45namespace LibUtilities
46{
47
48// =====================================================================
49// ==== DEFINITION OF THE CLASS NekSysOperators
50// ==== Defines operators needed by iterative solver of nonlinear or
51// ==== linear system
52// =====================================================================
54{
55public:
58
59 typedef std::function<void(InArrayType &, OutArrayType &, const bool &)>
61 typedef std::function<void(InArrayType &, InArrayType &, OutArrayType &,
62 const bool &)>
66 static const int nfunctor1 = 4;
67 static const int nfunctor2 = 1;
68
70 {
71 }
72
75 {
76 for (int i = 0; i < nfunctor1; ++i)
77 {
78 m_functors1[i] = in.m_functors1[i];
79 }
80 for (int i = 0; i < nfunctor2; ++i)
81 {
82 m_functors2[i] = in.m_functors2[i];
83 }
84 }
85
87 {
88 for (int i = 0; i < nfunctor1; ++i)
89 {
90 m_functors1[i] = in.m_functors1[i];
91 }
92 for (int i = 0; i < nfunctor2; ++i)
93 {
94 m_functors2[i] = in.m_functors2[i];
95 }
96
97 return *this;
98 }
99
100 template <typename FuncPointerT, typename ObjectPointerT>
101 void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
102 {
103 m_functors1[0] =
104 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
105 std::placeholders::_3);
106 }
107 template <typename FuncPointerT, typename ObjectPointerT>
108 void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
109 {
110 m_functors1[1] =
111 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
112 std::placeholders::_3);
113 }
114 template <typename FuncPointerT, typename ObjectPointerT>
115 void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
116 {
117 m_functors1[2] =
118 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
119 std::placeholders::_3);
120 }
121 template <typename FuncPointerT, typename ObjectPointerT>
122 void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
123 {
124 m_functors1[3] =
125 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
126 std::placeholders::_3);
127 }
128 template <typename FuncPointerT, typename ObjectPointerT>
129 void DefineNekSysFixPointIte(FuncPointerT func, ObjectPointerT obj)
130 {
131 m_functors2[0] =
132 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
133 std::placeholders::_3, std::placeholders::_4);
134 }
135
136 inline void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray,
137 const bool &flag = false) const
138 {
139 ASSERTL1(m_functors1[0], "DoNekSysResEval should be defined");
140 m_functors1[0](inarray, outarray, flag);
141 }
142
143 inline void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray,
144 const bool &flag = false) const
145 {
146 ASSERTL1(m_functors1[1], "DoNekSysLhsEval should be defined");
147 m_functors1[1](inarray, outarray, flag);
148 }
149
150 inline void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray,
151 const bool &flag = false) const
152 {
153 if (m_functors1[2])
154 {
155 m_functors1[2](inarray, outarray, flag);
156 }
157 else
158 {
159 Vmath::Vcopy(outarray.size(), inarray, 1, outarray, 1);
160 }
161 }
162
164 const bool &flag = false) const
165 {
166 ASSERTL1(m_functors1[3], "DoAssembleLoc should be defined");
167 m_functors1[3](xn, xn1, flag);
168 }
169
171 OutArrayType &xn1,
172 const bool &flag = false) const
173 {
174 ASSERTL1(m_functors2[0], "DoNekSysFixPointIte should be defined");
175 m_functors2[0](rhs, xn, xn1, flag);
176 }
177
178protected:
179 /* Defines operators
180 DoNekSysResEval :
181 evaluations the residual of the Nonlinear/Linear system
182 ie. the residual b-Ax and N(x) for linear and
183 nonlinear systems, respectively
184 May not be used for linear system.
185 DoNekSysLhsEval :
186 evaluations the LHS of the Nonlinear/Linear system (Ax),
187 where A is the matrix and x is solution vector.
188 For linear system A is the coefficient matrix;
189 For nonlinear system A is the coefficient matrix in
190 each nonlinear iterations, for example A is the
191 Jacobian matrix for Newton method;
192 DoNekSysPrecon :
193 Preconditioning operator of the system.
194 DoAssembleLoc :
195 Operator to assemble array and scater it back to the lcoal storage
196 if flag is true will zero Dirichlet conditions
197 DoNekSysFixPointIte :
198 */
200
201 /* Defines one operator
202 DoNekSysFixPointIte :
203 Operator to calculate RHS of fixed point iterations
204 (x^{n+1}=M^{-1}(b-N*x^{n}), with M+N=A).
205 */
207};
208
210{
211public:
213 {
214 }
215
217 {
218 }
219
231 bool m_DifferenceFlag0 = false;
232 bool m_DifferenceFlag1 = false;
233 bool m_useProjection = false;
234 std::string m_LinSysIterSolverTypeInNonlin = "GMRES";
235};
236
237class NekSys;
238
239typedef std::shared_ptr<NekSys> NekSysSharedPtr;
240
241class NekSys : public std::enable_shared_from_this<NekSys>
242{
243public:
244 /// Support creation through MemoryManager.
245 friend class MemoryManager<NekSys>;
246
249 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
250 const NekSysKey &pKey)
251 {
253 pSession, vRowComm, nDimen, pKey);
254 return p;
255 }
258 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
259 const NekSysKey &pKey);
261 {
262 v_InitObject();
263 }
265
267 {
268 m_operator = in;
269 }
270
272 {
273 return m_operator;
274 }
275
277 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
278 Array<OneD, NekDouble> &pOutput, const int nDir,
279 const NekDouble tol = 1.0E-7, const NekDouble factor = 1.0)
280 {
281 return v_SolveSystem(nGlobal, pInput, pOutput, nDir, tol, factor);
282 }
283
285 const int nIteration, const Array<OneD, const NekDouble> &Residual,
286 const NekDouble tol = 1.0E-7)
287 {
288 return v_ConvergenceCheck(nIteration, Residual, tol);
289 }
290
292 {
293 m_FlagWarnings = in;
294 }
295
296protected:
297 /// Maximum iterations
299 /// Tolerance of iterative solver.
301 /// Communicate.
303 /// Whether the iteration has been converged
305 /// Root if parallel
306 bool m_root;
307 /// Verbose
310 /// Operators
312 /// The dimension of the system
314
315 virtual void v_InitObject()
316 {
317 }
318
319 virtual int v_SolveSystem(const int nGlobal,
320 const Array<OneD, const NekDouble> &pInput,
321 Array<OneD, NekDouble> &pOutput, const int nDir,
322 const NekDouble tol, const NekDouble factor)
323 {
324 boost::ignore_unused(nGlobal, pInput, pOutput, nDir, tol, factor);
325 ASSERTL0(false, "LinSysIterSolver NOT CORRECT.");
326 return 0;
327 }
328
329 virtual bool v_ConvergenceCheck(
330 const int nIteration, const Array<OneD, const NekDouble> &Residual,
331 const NekDouble tol);
332
334 const Array<OneD, const NekDouble> &pInput,
335 Array<OneD, NekDouble> &pguess);
336};
337} // namespace LibUtilities
338} // namespace Nektar
339#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define LIB_UTILITIES_EXPORT
void SetFlagWarnings(bool in)
Definition: NekSys.h:291
bool m_root
Root if parallel.
Definition: NekSys.h:306
const NekSysOperators & GetSysOperators()
Definition: NekSys.h:271
static NekSysSharedPtr CreateInstance(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Definition: NekSys.h:247
virtual void v_NekSysInitialGuess(const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pguess)
Definition: NekSys.cpp:97
void SetSysOperators(const NekSysOperators &in)
Definition: NekSys.h:266
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:302
virtual void v_InitObject()
Definition: NekSys.h:315
bool ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol=1.0E-7)
Definition: NekSys.h:284
bool m_verbose
Verbose.
Definition: NekSys.h:308
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:300
NekSysOperators m_operator
Operators.
Definition: NekSys.h:311
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:313
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:319
NekSys(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Definition: NekSys.cpp:50
virtual bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol)
Definition: NekSys.cpp:76
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:276
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:304
int m_maxiter
Maximum iterations.
Definition: NekSys.h:298
NekDouble m_NonlinIterTolRelativeL2
Definition: NekSys.h:225
NekDouble m_LinSysRelativeTolInNonlin
Definition: NekSys.h:226
std::string m_LinSysIterSolverTypeInNonlin
Definition: NekSys.h:234
void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:122
void DefineNekSysFixPointIte(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:129
Array< OneD, FunctorType2 > FunctorType2Array
Definition: NekSys.h:65
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:136
void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:101
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:108
const Array< OneD, const NekDouble > InArrayType
Definition: NekSys.h:56
NekSysOperators & operator=(const NekSysOperators &in)
Definition: NekSys.h:86
std::function< void(InArrayType &, InArrayType &, OutArrayType &, const bool &)> FunctorType2
Definition: NekSys.h:63
Array< OneD, NekDouble > OutArrayType
Definition: NekSys.h:57
std::function< void(InArrayType &, OutArrayType &, const bool &)> FunctorType1
Definition: NekSys.h:60
NekSysOperators(const NekSysOperators &in)
Definition: NekSys.h:73
void DoNekSysFixPointIte(InArrayType &rhs, InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:170
void DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:163
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:150
Array< OneD, FunctorType1 > FunctorType1Array
Definition: NekSys.h:64
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:115
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:143
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:239
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
static const NekDouble kNekIterativeTol
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191