integrationScheme.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::integrationScheme
28 
29 Description
30  Base for a set of schemes which integrate simple ODEs which arise from
31  semi-implcit rate expressions.
32 
33  \f[
34  \frac{d \phi}{d t} = A - B \phi
35  \f]
36 
37  The methods are defined in terms of the effective time-step \f$\Delta
38  t_e\f$ by which the explicit rate is multiplied. The effective time-step is
39  a function of the actual time step and the implicit coefficient, which must
40  be implemented in each derived scheme.
41 
42  \f[
43  \Delta t_e = f(\Delta t, B)
44  \f]
45  \f[
46  \Delta \phi = (A - B \phi^n) \Delta t_e
47  \f]
48 
49  This class also facilitates integration in stages. If the explicit and
50  implicit coefficients, \f$A\f$ and \f$B\f$, are a summation of differing
51  contributions, \f$\sum \alpha_i\f$ and \f$\sum \beta_i\f$, then the
52  integration can be split up to determine the effect of each contribution.
53 
54  \f[
55  \frac{d \phi_i}{d t} = \alpha_i - \beta_i \phi
56  \f]
57  \f[
58  \Delta \phi_i = \alpha_i \Delta t -
59  \beta_i \int_0^{\Delta t} \phi d t
60  \f]
61  \f[
62  \Delta \phi_i = (\alpha_i - \beta_i \phi^n) \Delta t -
63  (A - B \phi^n) \int_0^{\Delta t} t_e dt
64  \f]
65 
66  These partial calculations are defined in terms of the integral of the
67  effective time-step, \f$\int_0^{\Delta t} t_e dt\f$, which is also
68  implemented in every derivation.
69 
70 SourceFiles
71  integrationScheme.C
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef integrationScheme_H
76 #define integrationScheme_H
77 
78 #include "autoPtr.H"
79 #include "runTimeSelectionTables.H"
80 #include "dictionary.H"
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 
87 /*---------------------------------------------------------------------------*\
88  Class integrationScheme Declaration
89 \*---------------------------------------------------------------------------*/
90 
92 {
93 public:
94 
95  //- Runtime type information
96  TypeName("integrationScheme");
97 
98 
99  //- Declare runtime constructor selection table
100 
102  (
103  autoPtr,
105  word,
106  (),
107  ()
108  );
109 
110 
111  // Constructors
112 
113  //- Construct
115 
116  //- Construct and return clone
117  virtual autoPtr<integrationScheme> clone() const = 0;
118 
119 
120  // Selectors
121 
122  //- Select an integration scheme
124  (
125  const word& phiName,
126  const dictionary& dict
127  );
128 
129 
130  //- Destructor
131  virtual ~integrationScheme();
132 
133 
134  // Member Functions
135 
136  //- Perform the integration explicitly
137  template<class Type>
138  static Type explicitDelta
139  (
140  const Type& phi,
141  const scalar dtEff,
142  const Type& Alpha,
143  const scalar Beta
144  );
145 
146  //- Perform the integration
147  template<class Type>
148  Type delta
149  (
150  const Type& phi,
151  const scalar dt,
152  const Type& Alpha,
153  const scalar Beta
154  ) const;
155 
156  //- Perform a part of the integration
157  template<class Type>
158  Type partialDelta
159  (
160  const Type& phi,
161  const scalar dt,
162  const Type& Alpha,
163  const scalar Beta,
164  const Type& alphai,
165  const scalar betai
166  ) const;
167 
168  //- Return the integration effective time step
169  virtual scalar dtEff(const scalar dt, const scalar Beta) const = 0;
170 
171  //- Return the integral of the effective time step
172  virtual scalar sumDtEff(const scalar dt, const scalar Beta) const = 0;
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace Foam
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 #ifdef NoRepository
184 #endif
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
Foam::integrationScheme::delta
Type delta(const Type &phi, const scalar dt, const Type &Alpha, const scalar Beta) const
Perform the integration.
Definition: integrationSchemeTemplates.C:47
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::integrationScheme::dtEff
virtual scalar dtEff(const scalar dt, const scalar Beta) const =0
Return the integration effective time step.
Foam::integrationScheme::integrationScheme
integrationScheme()
Construct.
Definition: integrationScheme.C:40
Foam::integrationScheme
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
Definition: integrationScheme.H:90
Foam::integrationScheme::explicitDelta
static Type explicitDelta(const Type &phi, const scalar dtEff, const Type &Alpha, const scalar Beta)
Perform the integration explicitly.
Definition: integrationSchemeTemplates.C:34
Foam::integrationScheme::sumDtEff
virtual scalar sumDtEff(const scalar dt, const scalar Beta) const =0
Return the integral of the effective time step.
Foam::integrationScheme::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, integrationScheme, word,(),())
Declare runtime constructor selection table.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::integrationScheme::~integrationScheme
virtual ~integrationScheme()
Destructor.
Definition: integrationScheme.C:46
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
integrationSchemeTemplates.C
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::integrationScheme::New
static autoPtr< integrationScheme > New(const word &phiName, const dictionary &dict)
Select an integration scheme.
Definition: integrationSchemeNew.C:35
Foam::integrationScheme::partialDelta
Type partialDelta(const Type &phi, const scalar dt, const Type &Alpha, const scalar Beta, const Type &alphai, const scalar betai) const
Perform a part of the integration.
Definition: integrationSchemeTemplates.C:60
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::integrationScheme::TypeName
TypeName("integrationScheme")
Runtime type information.
dictionary.H
Foam::integrationScheme::clone
virtual autoPtr< integrationScheme > clone() const =0
Construct and return clone.
autoPtr.H