ODESolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "ODESolver.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(ODESolver, 0);
36  defineRunTimeSelectionTable(ODESolver, dictionary);
37 }
38 
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 (
44  const scalarField& y0,
45  const scalarField& y,
46  const scalarField& err
47 ) const
48 {
49  // Calculate the maximum error
50  scalar maxErr = 0.0;
51  forAll(err, i)
52  {
53  scalar tol = absTol_[i] + relTol_[i]*max(mag(y0[i]), mag(y[i]));
54  maxErr = max(maxErr, mag(err[i])/tol);
55  }
56 
57  return maxErr;
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 :
65  odes_(ode),
66  maxN_(ode.nEqns()),
67  n_(ode.nEqns()),
68  absTol_(n_, dict.getOrDefault<scalar>("absTol", SMALL)),
69  relTol_(n_, dict.getOrDefault<scalar>("relTol", 1e-4)),
70  maxSteps_(dict.getOrDefault<label>("maxSteps", 10000))
71 {}
72 
73 
75 (
76  const ODESystem& ode,
77  const scalarField& absTol,
78  const scalarField& relTol
79 )
80 :
81  odes_(ode),
82  maxN_(ode.nEqns()),
83  n_(ode.nEqns()),
84  absTol_(absTol),
85  relTol_(relTol),
86  maxSteps_(10000)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  if (odes_.nEqns() != n_)
95  {
96  if (odes_.nEqns() > maxN_)
97  {
99  << "Specified number of equations " << odes_.nEqns()
100  << " greater than maximum " << maxN_
101  << abort(FatalError);
102  }
103 
104  n_ = odes_.nEqns();
105 
106  resizeField(absTol_);
107  resizeField(relTol_);
108 
109  return true;
110  }
111 
112  return false;
113 }
114 
115 
117 (
118  scalar& x,
119  scalarField& y,
120  scalar& dxTry
121 ) const
122 {
123  stepState step(dxTry);
124  solve(x, y, step);
125  dxTry = step.dxTry;
126 }
127 
128 
130 (
131  scalar& x,
132  scalarField& y,
133  stepState& step
134 ) const
135 {
136  scalar x0 = x;
137  solve(x, y, step.dxTry);
138  step.dxDid = x - x0;
139 }
140 
141 
143 (
144  const scalar xStart,
145  const scalar xEnd,
146  scalarField& y,
147  scalar& dxTry
148 ) const
149 {
150  stepState step(dxTry);
151  scalar x = xStart;
152 
153  for (label nStep=0; nStep<maxSteps_; ++nStep)
154  {
155  // Store previous iteration dxTry
156  scalar dxTry0 = step.dxTry;
157 
158  step.reject = false;
159 
160  // Check if this is a truncated step and set dxTry to integrate to xEnd
161  if ((x + step.dxTry - xEnd)*(x + step.dxTry - xStart) > 0)
162  {
163  step.last = true;
164  step.dxTry = xEnd - x;
165  }
166 
167  // Integrate as far as possible up to step.dxTry
168  solve(x, y, step);
169 
170  // Check if reached xEnd
171  if ((x - xEnd)*(xEnd - xStart) >= 0)
172  {
173  if (nStep > 0 && step.last)
174  {
175  step.dxTry = dxTry0;
176  }
177 
178  dxTry = step.dxTry;
179 
180  return;
181  }
182 
183  step.first = false;
184 
185  // If the step.dxTry was reject set step.prevReject
186  if (step.reject)
187  {
188  step.prevReject = true;
189  }
190  }
191 
193  << "Integration steps greater than maximum " << maxSteps_ << nl
194  << " xStart = " << xStart << ", xEnd = " << xEnd
195  << ", x = " << x << ", dxDid = " << step.dxDid << nl
196  << " y = " << y
197  << exit(FatalError);
198 }
199 
200 
201 // ************************************************************************* //
Foam::ODESolver::normalizeError
scalar normalizeError(const scalarField &y0, const scalarField &y, const scalarField &err) const
Return the nomalized scalar error.
Definition: ODESolver.C:43
Foam::ODESolver::stepState
Definition: ODESolver.H:106
Foam::ODESolver::stepState::reject
bool reject
Definition: ODESolver.H:115
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field< scalar >
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:281
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::ode
An ODE solver for chemistry.
Definition: ode.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::ODESolver::stepState::dxTry
scalar dxTry
Definition: ODESolver.H:111
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
ODESolver.H
Foam::ODESolver::stepState::first
bool first
Definition: ODESolver.H:113
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::ODESolver::ODESolver
ODESolver(const ODESolver &)=delete
No copy construct.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::ODESolver::stepState::last
bool last
Definition: ODESolver.H:114
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::ODESystem
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:49
Foam::ODESolver::solve
virtual void solve(scalar &x, scalarField &y, scalar &dxTry) const
Solve the ODE system as far as possible up to dxTry.
Definition: ODESolver.C:117
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::ODESolver::stepState::dxDid
scalar dxDid
Definition: ODESolver.H:112
Foam::ODESolver::resize
virtual bool resize()=0
Resize the ODE solver.
Definition: ODESolver.C:92
Foam::ODESolver::stepState::prevReject
bool prevReject
Definition: ODESolver.H:116
y
scalar y
Definition: LISASMDCalcMethod1.H:14