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