CrankNicolsonDdtScheme.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::fv::CrankNicolsonDdtScheme
29
30Group
31 grpFvDdtSchemes
32
33Description
34 Second-oder Crank-Nicolson implicit ddt using the current and
35 previous time-step fields as well as the previous time-step ddt.
36
37 The Crank-Nicolson scheme is often unstable for complex flows in complex
38 geometries and it is necessary to "off-centre" the scheme to stabilize it
39 while retaining greater temporal accuracy than the first-order
40 Euler-implicit scheme. Off-centering is specified via the mandatory
41 coefficient \c ocCoeff in the range [0,1] following the scheme name e.g.
42 \verbatim
43 ddtSchemes
44 {
45 default CrankNicolson 0.9;
46 }
47 \endverbatim
48 or with an optional "ramp" function to transition from the Euler scheme to
49 Crank-Nicolson over a initial period to avoid start-up problems, e.g.
50 \verbatim
51 ddtSchemes
52 {
53 default CrankNicolson
54 ocCoeff
55 {
56 type scale;
57 scale linearRamp;
58 duration 0.01;
59 value 0.9;
60 };
61 }
62 \endverbatim
63 With a coefficient of 1 the scheme is fully centred and second-order,
64 with a coefficient of 0 the scheme is equivalent to Euler-implicit.
65 A coefficient of 0.9 has been found to be suitable for a range of cases for
66 which higher-order accuracy in time is needed and provides similar accuracy
67 and stability to the backward scheme. However, the backward scheme has
68 been found to be more robust and provides formal second-order accuracy in
69 time.
70
71 The advantage of the Crank-Nicolson scheme over backward is that only the
72 new and old-time values are needed, the additional terms relating to the
73 fluxes and sources are evaluated at the mid-point of the time-step which
74 provides the opportunity to limit the fluxes in such a way as to ensure
75 boundedness while maintaining greater accuracy in time compared to the
76 Euler-implicit scheme. This approach is now used with MULES in the
77 interFoam family of solvers. Boundedness cannot be guaranteed with the
78 backward scheme.
79
80Note
81 The Crank-Nicolson coefficient for the implicit part of the RHS is related
82 to the off-centering coefficient by
83 \verbatim
84 cnCoeff = 1.0/(1.0 + ocCoeff);
85 \endverbatim
86
87See also
88 Foam::Euler
89 Foam::backward
90
91SourceFiles
92 CrankNicolsonDdtScheme.C
93
94\*---------------------------------------------------------------------------*/
95
96#ifndef CrankNicolsonDdtScheme_H
97#define CrankNicolsonDdtScheme_H
98
99#include "ddtScheme.H"
100#include "Function1.H"
101
102// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103
104namespace Foam
105{
106
107// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108
109namespace fv
110{
111
112/*---------------------------------------------------------------------------*\
113 Class CrankNicolsonDdtScheme Declaration
114\*---------------------------------------------------------------------------*/
115
116template<class Type>
118:
119 public fv::ddtScheme<Type>
120{
121 // Private Data
122
123 //- Class to store the ddt0 fields on the objectRegistry for use in the
124 // next time-step. The start-time index of the CN scheme is also
125 // stored to help handle the transition from Euler to CN
126 template<class GeoField>
127 class DDt0Field
128 :
129 public GeoField
130 {
131 label startTimeIndex_;
132
133 public:
134
135 //- Constructor from file for restart.
136 DDt0Field
137 (
138 const IOobject& io,
139 const fvMesh& mesh
140 );
141
142 //- Constructor from components, initisalised to zero with given
143 // dimensions.
144 DDt0Field
145 (
146 const IOobject& io,
147 const fvMesh& mesh,
149 );
150
151 //- Return the start-time index
152 label startTimeIndex() const;
153
154 //- Cast to the underlying GeoField
155 GeoField& operator()();
156
157 //- Assignment to a GeoField
158 void operator=(const GeoField& gf);
159 };
160
161
162 //- Off-centering coefficient function
163 // 1 -> CN, less than one blends with EI
165
166
167 // Private Member Functions
168
169 //- No copy construct
171
172 //- No copy assignment
173 void operator=(const CrankNicolsonDdtScheme&) = delete;
174
175 template<class GeoField>
176 DDt0Field<GeoField>& ddt0_
177 (
178 const word& name,
179 const dimensionSet& dims
180 );
181
182 //- Check if the ddt0 needs to be evaluated for this time-step
183 template<class GeoField>
184 bool evaluate(DDt0Field<GeoField>& ddt0);
185
186 //- Return the coefficient for Euler scheme for the first time-step
187 // for and CN thereafter
188 template<class GeoField>
189 scalar coef_(const DDt0Field<GeoField>&) const;
190
191 //- Return the old time-step coefficient for Euler scheme for the
192 // second time-step and for CN thereafter
193 template<class GeoField>
194 scalar coef0_(const DDt0Field<GeoField>&) const;
195
196 //- Return the reciprocal time-step coefficient for Euler for the
197 // first time-step and CN thereafter
198 template<class GeoField>
199 dimensionedScalar rDtCoef_(const DDt0Field<GeoField>&) const;
200
201 //- Return the reciprocal old time-step coefficient for Euler for the
202 // second time-step and CN thereafter
203 template<class GeoField>
204 dimensionedScalar rDtCoef0_(const DDt0Field<GeoField>&) const;
205
206 //- Return ddt0 multiplied by the off-centreing coefficient
207 template<class GeoField>
208 tmp<GeoField> offCentre_(const GeoField& ddt0) const;
209
210
211public:
212
213 //- Runtime type information
214 TypeName("CrankNicolson");
215
216
217 // Constructors
218
219 //- Construct from mesh
221
222 //- Construct from mesh and Istream
224
225
226 // Member Functions
227
228 //- Return mesh reference
229 const fvMesh& mesh() const
230 {
232 }
233
234 //- Return the current off-centreing coefficient
235 scalar ocCoeff() const
236 {
237 return ocCoeff_->value(mesh().time().value());
238 }
239
241 (
242 const dimensioned<Type>&
243 );
244
246 (
248 );
249
251 (
252 const dimensionedScalar&,
254 );
255
257 (
258 const volScalarField&,
260 );
261
263 (
264 const volScalarField& alpha,
265 const volScalarField& rho,
267 );
268
270 (
272 );
273
275 (
276 const dimensionedScalar&,
278 );
279
281 (
282 const volScalarField&,
284 );
285
287 (
288 const volScalarField& alpha,
289 const volScalarField& rho,
291 );
294
296 (
299 );
300
302 (
304 const fluxFieldType& phi
305 );
306
308 (
309 const volScalarField& rho,
312 );
313
315 (
316 const volScalarField& rho,
318 const fluxFieldType& phi
319 );
320
321
323 (
325 );
326};
327
328
329template<>
331(
334);
335
336template<>
338(
339 const volScalarField& U,
341);
342
343template<>
345(
346 const volScalarField& rho,
347 const volScalarField& U,
349);
350
351template<>
353(
354 const volScalarField& rho,
355 const volScalarField& U,
357);
358
359
360// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361
362} // End namespace fv
363
364// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365
366} // End namespace Foam
367
368// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369
370#ifdef NoRepository
371 #include "CrankNicolsonDdtScheme.C"
372#endif
373
374// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375
376#endif
377
378// ************************************************************************* //
surfaceScalarField & phi
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
scalar ocCoeff() const
Return the current off-centreing coefficient.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
TypeName("CrankNicolson")
Runtime type information.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:87
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const volScalarField & psi
autoPtr< surfaceVectorField > Uf
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
labelList fv(nPoints)
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73