ODESolver.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::ODESolver
26 
27 Description
28  Abstract base-class for ODE system solvers
29 
30 SourceFiles
31  ODESolver.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef ODESolver_H
36 #define ODESolver_H
37 
38 #include "ODESystem.H"
39 #include "typeInfo.H"
40 #include "autoPtr.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class ODESolver Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class ODESolver
52 {
53 
54 protected:
55 
56  // Protected data
57 
58  //- Reference to ODESystem
59  const ODESystem& odes_;
60 
61  //- Maximum size of the ODESystem
62  const label maxN_;
63 
64  //- Size of the ODESystem (adjustable)
65  mutable label n_;
66 
67  //- Absolute convergence tolerance per step
69 
70  //- Relative convergence tolerance per step
72 
73  //- The maximum number of sub-steps allowed for the integration step
75 
76 
77  // Protected Member Functions
78 
79  //- Return the normalised scalar error
80  scalar normaliseError
81  (
82  const scalarField& y0,
83  const scalarField& y,
84  const scalarField& err
85  ) const;
86 
87 
88 public:
89 
90  friend class ODESystem;
91 
92  //- Runtime type information
93  TypeName("ODESolver");
94 
95  class stepState
96  {
97  public:
98 
99  const bool forward;
100  scalar dxTry;
101  scalar dxDid;
102  bool first;
103  bool last;
104  bool reject;
105  bool prevReject;
107  stepState(const scalar dx)
108  :
109  forward(dx > 0 ? true : false),
110  dxTry(dx),
111  dxDid(0),
112  first(true),
113  last(false),
114  reject(false),
115  prevReject(false)
116  {}
117  };
118 
119 
120  // Declare run-time constructor selection table
121 
123  (
124  autoPtr,
125  ODESolver,
126  dictionary,
127  (const ODESystem& ode, const dictionary& dict),
128  (ode, dict)
129  );
130 
131 
132  // Constructors
133 
134  //- Construct for given ODESystem
135  ODESolver(const ODESystem& ode, const dictionary& dict);
136 
137  //- Construct for given ODESystem specifying tolerances
138  ODESolver
139  (
140  const ODESystem& ode,
141  const scalarField& absTol,
142  const scalarField& relTol
143  );
144 
145  //- Disallow default bitwise copy construction
146  ODESolver(const ODESolver&) = delete;
147 
148 
149  // Selectors
150 
151  //- Select null constructed
152  static autoPtr<ODESolver> New
153  (
154  const ODESystem& ode,
155  const dictionary& dict
156  );
157 
158 
159  //- Destructor
160  virtual ~ODESolver()
161  {}
162 
163 
164  // Member Functions
165 
166  //- Return the number of equations to solve
167  inline label nEqns() const;
168 
169  //- Return access to the absolute tolerance field
170  inline scalarField& absTol();
171 
172  //- Return access to the relative tolerance field
173  inline scalarField& relTol();
174 
175  //- Resize the ODE solver
176  virtual bool resize() = 0;
177 
178  template<class Type>
179  static inline void resizeField(UList<Type>& f, const label n);
180 
181  template<class Type>
182  inline void resizeField(UList<Type>& f) const;
183 
184  inline void resizeMatrix(scalarSquareMatrix& m) const;
185 
186  //- Solve the ODE system from the current state xStart, y
187  // and the optional index into the list of systems to solve li
188  // as far as possible up to dxTry adjusting the step as necessary
189  // to provide a solution within the specified tolerance.
190  // Update the state and return an estimate for the next step in dxTry
191  virtual void solve
192  (
193  scalar& x,
194  scalarField& y,
195  const label li,
196  scalar& dxTry
197  ) const;
198 
199  //- Solve the ODE system from the current state xStart, y
200  // and the optional index into the list of systems to solve li
201  // as far as possible up to step.dxTry adjusting the step as necessary
202  // to provide a solution within the specified tolerance.
203  // Update the state and return the step actually taken in step.dxDid
204  // and an estimate for the next step in step.dxTry
205  virtual void solve
206  (
207  scalar& x,
208  scalarField& y,
209  const label li,
210  stepState& step
211  ) const;
212 
213  //- Solve the ODE system from the current state xStart, y
214  // and the optional index into the list of systems to solve li
215  // to xEnd and return an estimate for the next step in dxTry
216  virtual void solve
217  (
218  const scalar xStart,
219  const scalar xEnd,
220  scalarField& y,
221  const label li,
222  scalar& dxEst
223  ) const;
224 
225 
226  // Member Operators
227 
228  //- Disallow default bitwise assignment
229  void operator=(const ODESolver&) = delete;
230 };
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 } // End namespace Foam
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #include "ODESolverI.H"
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #endif
244 
245 // ************************************************************************* //
label nEqns() const
Return the number of equations to solve.
Definition: ODESolverI.H:29
ODESolver(const ODESystem &ode, const dictionary &dict)
Construct for given ODESystem.
Definition: ODESolver.C:60
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
scalarField & relTol()
Return access to the relative tolerance field.
Definition: ODESolverI.H:41
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:46
virtual bool resize()=0
Resize the ODE solver.
Definition: ODESolver.C:89
void operator=(const ODESolver &)=delete
Disallow default bitwise assignment.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedScalar y0(const dimensionedScalar &ds)
An ODE solver for chemistry.
Definition: ode.H:50
stepState(const scalar dx)
Definition: ODESolver.H:106
scalar normaliseError(const scalarField &y0, const scalarField &y, const scalarField &err) const
Return the normalised scalar error.
Definition: ODESolver.C:40
virtual ~ODESolver()
Destructor.
Definition: ODESolver.H:159
scalar y
void resizeMatrix(scalarSquareMatrix &m) const
Definition: ODESolverI.H:61
const ODESystem & odes_
Reference to ODESystem.
Definition: ODESolver.H:58
TypeName("ODESolver")
Runtime type information.
label maxSteps_
The maximum number of sub-steps allowed for the integration step.
Definition: ODESolver.H:73
const label maxN_
Maximum size of the ODESystem.
Definition: ODESolver.H:61
labelList f(nPoints)
declareRunTimeSelectionTable(autoPtr, ODESolver, dictionary,(const ODESystem &ode, const dictionary &dict),(ode, dict))
virtual void solve(scalar &x, scalarField &y, const label li, scalar &dxTry) const
Solve the ODE system from the current state xStart, y.
Definition: ODESolver.C:116
scalarField relTol_
Relative convergence tolerance per step.
Definition: ODESolver.H:70
Abstract base-class for ODE system solvers.
Definition: ODESolver.H:50
scalarField absTol_
Absolute convergence tolerance per step.
Definition: ODESolver.H:67
static autoPtr< ODESolver > New(const ODESystem &ode, const dictionary &dict)
Select null constructed.
Definition: ODESolverNew.C:31
scalarField & absTol()
Return access to the absolute tolerance field.
Definition: ODESolverI.H:35
label n
static void resizeField(UList< Type > &f, const label n)
Definition: ODESolverI.H:48
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
label n_
Size of the ODESystem (adjustable)
Definition: ODESolver.H:64
Namespace for OpenFOAM.