ODESolver.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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  //- Size of the ODESystem
62  label n_;
63 
64  //- Absolute convergence tolerance per step
66 
67  //- Relative convergence tolerance per step
69 
70  //- The maximum number of sub-steps allowed for the integration step
72 
73 
74  // Protected Member Functions
75 
76  //- Return the nomalized scalar error
77  scalar normalizeError
78  (
79  const scalarField& y0,
80  const scalarField& y,
81  const scalarField& err
82  ) const;
83 
84  //- Disallow default bitwise copy construct
85  ODESolver(const ODESolver&);
86 
87  //- Disallow default bitwise assignment
88  void operator=(const ODESolver&);
89 
90 
91 public:
92 
93  //- Runtime type information
94  TypeName("ODESolver");
95 
96  class stepState
97  {
98  public:
99 
100  const bool forward;
101  scalar dxTry;
102  scalar dxDid;
103  bool first;
104  bool last;
105  bool reject;
106  bool prevReject;
108  stepState(const scalar dx)
109  :
110  forward(dx > 0 ? true : false),
111  dxTry(dx),
112  dxDid(0),
113  first(true),
114  last(false),
115  reject(false),
116  prevReject(false)
117  {}
118  };
119 
120 
121  // Declare run-time constructor selection table
122 
124  (
125  autoPtr,
126  ODESolver,
127  dictionary,
128  (const ODESystem& ode, const dictionary& dict),
129  (ode, dict)
130  );
131 
132 
133  // Constructors
134 
135  //- Construct for given ODESystem
136  ODESolver(const ODESystem& ode, const dictionary& dict);
137 
138  //- Construct for given ODESystem specifying tolerances
139  ODESolver
140  (
141  const ODESystem& ode,
142  const scalarField& absTol,
143  const scalarField& relTol
144  );
145 
146 
147  // Selectors
148 
149  //- Select null constructed
150  static autoPtr<ODESolver> New
151  (
152  const ODESystem& ode,
153  const dictionary& dict
154  );
155 
156 
157  //- Destructor
158  virtual ~ODESolver()
159  {}
160 
161 
162  // Member Functions
165  {
166  return absTol_;
167  }
170  {
171  return relTol_;
172  }
173 
174  //- Solve the ODE system as far as possible upto dxTry
175  // adjusting the step as necessary to provide a solution within
176  // the specified tolerance.
177  // Update the state and return an estimate for the next step in dxTry
178  virtual void solve
179  (
180  scalar& x,
181  scalarField& y,
182  scalar& dxTry
183  ) const //= 0;
184  {
185  stepState step(dxTry);
186  solve(x, y, step);
187  dxTry = step.dxTry;
188  }
189 
190  //- Solve the ODE system as far as possible upto dxTry
191  // adjusting the step as necessary to provide a solution within
192  // the specified tolerance.
193  // Update the state and return an estimate for the next step in dxTry
194  virtual void solve
195  (
196  scalar& x,
197  scalarField& y,
198  stepState& step
199  ) const;
200 
201  //- Solve the ODE system from xStart to xEnd, update the state
202  // and return an estimate for the next step in dxTry
203  virtual void solve
204  (
205  const scalar xStart,
206  const scalar xEnd,
207  scalarField& y,
208  scalar& dxEst
209  ) const;
210 };
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #endif
220 
221 // ************************************************************************* //
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
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:46
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedScalar y0(const dimensionedScalar &ds)
An ODE solver for chemistry.
Definition: ode.H:50
stepState(const scalar dx)
Definition: ODESolver.H:107
virtual ~ODESolver()
Destructor.
Definition: ODESolver.H:157
scalar y
scalarField & absTol()
Definition: ODESolver.H:163
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:70
ODESolver(const ODESolver &)
Disallow default bitwise copy construct.
virtual void solve(scalar &x, scalarField &y, scalar &dxTry) const
Solve the ODE system as far as possible upto dxTry.
Definition: ODESolver.H:178
declareRunTimeSelectionTable(autoPtr, ODESolver, dictionary,(const ODESystem &ode, const dictionary &dict),(ode, dict))
scalarField relTol_
Relative convergence tolerance per step.
Definition: ODESolver.H:67
scalar normalizeError(const scalarField &y0, const scalarField &y, const scalarField &err) const
Return the nomalized scalar error.
Definition: ODESolver.C:67
Abstract base-class for ODE system solvers.
Definition: ODESolver.H:50
void operator=(const ODESolver &)
Disallow default bitwise assignment.
scalarField absTol_
Absolute convergence tolerance per step.
Definition: ODESolver.H:64
static autoPtr< ODESolver > New(const ODESystem &ode, const dictionary &dict)
Select null constructed.
Definition: ODESolverNew.C:31
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
scalarField & relTol()
Definition: ODESolver.H:168
label n_
Size of the ODESystem.
Definition: ODESolver.H:61
Namespace for OpenFOAM.