localEulerDdtScheme.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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::fv::localEulerDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 Local time-step first-order Euler implicit/explicit ddt.
34
35 The reciprocal of the local time-step field is looked-up from the database.
36
37 This scheme should only be used for steady-state computations using
38 transient codes where local time-stepping is preferably to under-relaxation
39 for transport consistency reasons.
40
41See also
42 Foam::fv::CoEulerDdtScheme
43
44SourceFiles
45 localEulerDdt.C
46 localEulerDdtScheme.C
47 localEulerDdtSchemes.C
48
49\*---------------------------------------------------------------------------*/
50
51#ifndef localEulerDdtScheme_H
52#define localEulerDdtScheme_H
53
54#include "ddtScheme.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace fv
64{
65
66/*---------------------------------------------------------------------------*\
67 Class localEulerDdt Declaration
68\*---------------------------------------------------------------------------*/
70class localEulerDdt
71{
72public:
73
74 //- Name of the reciprocal local time-step field
75 static word rDeltaTName;
76
77 //- Name of the reciprocal local face time-step field
78 static word rDeltaTfName;
79
80 //- Name of the reciprocal local sub-cycling time-step field
81 static word rSubDeltaTName;
82
83
84 // Constructors
87 {}
88
89
90 // Member Functions
91
92 //- Return true if LTS is enabled
93 static bool enabled(const fvMesh& mesh);
94
95 //- Return the reciprocal of the local time-step
96 // looked-up from the objectRegistry
97 static const volScalarField& localRDeltaT(const fvMesh& mesh);
98
99 //- Return the reciprocal of the local face time-step
100 // looked-up from the objectRegistry
101 static const surfaceScalarField& localRDeltaTf(const fvMesh& mesh);
102
103 //- Calculate and return the reciprocal of the local sub-cycling
104 // time-step
106 (
107 const fvMesh& mesh,
108 const label nAlphaSubCycles
109 );
110};
111
112
113/*---------------------------------------------------------------------------*\
114 Class localEulerDdtScheme Declaration
115\*---------------------------------------------------------------------------*/
116
117template<class Type>
119:
120 public localEulerDdt,
121 public fv::ddtScheme<Type>
122{
123 // Private Member Functions
124
125 //- No copy construct
127
128 //- No copy assignment
129 void operator=(const localEulerDdtScheme&) = delete;
130
131 //- Return the reciprocal of the local time-step
132 const volScalarField& localRDeltaT() const;
133
134 //- Return the reciprocal of the local face time-step
135 const surfaceScalarField& localRDeltaTf() const;
136
137
138public:
139
140 //- Runtime type information
141 TypeName("localEuler");
142
143
144 // Constructors
145
146 //- Construct from mesh
148 :
149 ddtScheme<Type>(mesh)
150 {}
151
152 //- Construct from mesh and Istream
154 :
155 ddtScheme<Type>(mesh, is)
156 {}
157
158
159 // Member Functions
160
161 //- Return mesh reference
162 const fvMesh& mesh() const
163 {
165 }
166
168 (
169 const dimensioned<Type>&
170 );
171
173 (
175 );
176
178 (
179 const dimensionedScalar&,
181 );
182
184 (
185 const volScalarField&,
187 );
188
190 (
191 const volScalarField& alpha,
192 const volScalarField& rho,
194 );
195
197 (
199 );
200
202 (
204 );
205
207 (
208 const dimensionedScalar&,
210 );
211
213 (
214 const volScalarField&,
216 );
217
219 (
220 const volScalarField& alpha,
221 const volScalarField& rho,
223 );
226
228 (
231 );
232
234 (
236 const fluxFieldType& phi
237 );
238
240 (
241 const volScalarField& rho,
244 );
245
247 (
248 const volScalarField& rho,
250 const fluxFieldType& phi
251 );
252
254 (
256 );
257};
258
259
260template<>
262(
265);
266
267template<>
269(
270 const volScalarField& U,
272);
273
274template<>
276(
277 const volScalarField& rho,
278 const volScalarField& U,
280);
281
282template<>
284(
285 const volScalarField& rho,
286 const volScalarField& U,
288);
289
290
291// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292
293} // End namespace fv
294
295// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296
297} // End namespace Foam
298
299// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300
301#ifdef NoRepository
302 #include "localEulerDdtScheme.C"
303#endif
304
305// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307#endif
308
309// ************************************************************************* //
surfaceScalarField & phi
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for ddt schemes.
Definition: ddtScheme.H:87
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
Local time-step first-order Euler implicit/explicit ddt.
localEulerDdtScheme(const fvMesh &mesh)
Construct from mesh.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
TypeName("localEuler")
Runtime type information.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
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 > &)
localEulerDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
static word rSubDeltaTName
Name of the reciprocal local sub-cycling time-step field.
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:60
static word rDeltaTfName
Name of the reciprocal local face time-step field.
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:72
static word rDeltaTName
Name of the reciprocal local time-step field.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:48
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:39
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
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
Namespace for OpenFOAM.
labelList fv(nPoints)
volScalarField & alpha
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73