CoEulerDdtScheme.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::CoEulerDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 Courant number limited first-order Euler implicit/explicit ddt.
34
35 The time-step is adjusted locally so that the local Courant number
36 does not exceed the specified value.
37
38 This scheme should only be used for steady-state computations
39 using transient codes where local time-stepping is preferable to
40 under-relaxation for transport consistency reasons.
41
42SourceFiles
43 CoEulerDdtScheme.C
44
45\*---------------------------------------------------------------------------*/
46
47#ifndef CoEulerDdtScheme_H
48#define CoEulerDdtScheme_H
49
50#include "ddtScheme.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace fv
60{
61
62/*---------------------------------------------------------------------------*\
63 Class CoEulerDdtScheme Declaration
64\*---------------------------------------------------------------------------*/
65
66template<class Type>
68:
69 public fv::ddtScheme<Type>
70{
71 // Private Data
72
73 //- Name of the flux field used to calculate the local time-step
74 word phiName_;
75
76 //- Name of the density field used to obtain the volumetric flux
77 // from the mass flux if required
78 word rhoName_;
79
80 //- Maximum local Courant number
81 scalar maxCo_;
82
83
84 // Private Member Functions
85
86 //- No copy construct
87 CoEulerDdtScheme(const CoEulerDdtScheme&) = delete;
88
89 //- No copy assignment
90 void operator=(const CoEulerDdtScheme&) = delete;
91
92 //- Return the reciprocal of the Courant-number limited time-step
93 tmp<volScalarField> CorDeltaT() const;
94
95 //- Return the reciprocal of the face-Courant-number limited time-step
96 tmp<surfaceScalarField> CofrDeltaT() const;
97
98
99public:
100
101 //- Runtime type information
102 TypeName("CoEuler");
103
104
105 // Constructors
106
107 //- Construct from mesh and Istream
109 :
110 ddtScheme<Type>(mesh, is),
111 phiName_(is),
112 rhoName_(is),
113 maxCo_(readScalar(is))
114 {}
115
116
117 // Member Functions
118
119 //- Return mesh reference
120 const fvMesh& mesh() const
121 {
123 }
124
126 (
127 const dimensioned<Type>&
128 );
129
131 (
133 );
134
136 (
137 const dimensionedScalar&,
139 );
140
142 (
143 const volScalarField&,
145 );
146
148 (
149 const volScalarField& alpha,
150 const volScalarField& rho,
152 );
153
155 (
157 );
158
160 (
161 const dimensionedScalar&,
163 );
164
166 (
167 const volScalarField&,
169 );
170
172 (
173 const volScalarField& alpha,
174 const volScalarField& rho,
176 );
179
181 (
184 );
185
187 (
189 const fluxFieldType& phi
190 );
191
193 (
194 const volScalarField& rho,
197 );
198
200 (
201 const volScalarField& rho,
203 const fluxFieldType& phi
204 );
205
207 (
209 );
210};
211
212
213template<>
215(
218);
219
220template<>
222(
223 const volScalarField& U,
225);
226
227template<>
229(
230 const volScalarField& rho,
231 const volScalarField& U,
233);
234
235template<>
237(
238 const volScalarField& rho,
239 const volScalarField& U,
241);
242
243
244// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245
246} // End namespace fv
247
248// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249
250} // End namespace Foam
251
252// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253
254#ifdef NoRepository
255 #include "CoEulerDdtScheme.C"
256#endif
257
258// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259
260#endif
261
262// ************************************************************************* //
surfaceScalarField & phi
Generic GeometricField class.
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
Courant number limited first-order Euler implicit/explicit ddt.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
TypeName("CoEuler")
Runtime type information.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
CoEulerDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
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
Namespace for OpenFOAM.
labelList fv(nPoints)
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73