EulerDdtScheme.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::EulerDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 Basic first-order Euler implicit/explicit ddt using only the current and
34 previous time-step values.
35
36SourceFiles
37 EulerDdtScheme.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef EulerDdtScheme_H
42#define EulerDdtScheme_H
43
44#include "ddtScheme.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace fv
54{
55
56/*---------------------------------------------------------------------------*\
57 Class EulerDdtScheme Declaration
58\*---------------------------------------------------------------------------*/
59
60template<class Type>
62:
63 public ddtScheme<Type>
64{
65 // Private Member Functions
66
67 //- No copy construct
68 EulerDdtScheme(const EulerDdtScheme&) = delete;
69
70 //- No copy assignment
71 void operator=(const EulerDdtScheme&) = delete;
72
73
74public:
75
76 //- Runtime type information
77 TypeName("Euler");
78
79
80 // Constructors
81
82 //- Construct from mesh
84 :
85 ddtScheme<Type>(mesh)
86 {}
87
88 //- Construct from mesh and Istream
90 :
91 ddtScheme<Type>(mesh, is)
92 {}
93
94
95 // Member Functions
96
97 //- Return mesh reference
98 const fvMesh& mesh() const
99 {
101 }
102
104 (
105 const dimensioned<Type>&
106 );
107
109 (
111 );
112
114 (
115 const dimensionedScalar&,
117 );
118
120 (
121 const volScalarField&,
123 );
124
126 (
127 const volScalarField& alpha,
128 const volScalarField& rho,
130 );
131
133 (
135 );
136
138 (
140 );
141
143 (
144 const dimensionedScalar&,
146 );
147
149 (
150 const volScalarField&,
152 );
153
155 (
156 const volScalarField& alpha,
157 const volScalarField& rho,
159 );
162
164 (
167 );
168
170 (
172 const fluxFieldType& phi
173 );
174
176 (
177 const volScalarField& rho,
180 );
181
183 (
184 const volScalarField& rho,
186 const fluxFieldType& phi
187 );
188
190 (
192 );
193};
194
195
196template<>
198(
201);
202
203template<>
205(
206 const volScalarField& U,
208);
209
210template<>
212(
213 const volScalarField& rho,
214 const volScalarField& U,
216);
217
218template<>
220(
221 const volScalarField& rho,
222 const volScalarField& U,
224);
225
226
227// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228
229} // End namespace fv
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace Foam
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237#ifdef NoRepository
238 #include "EulerDdtScheme.C"
239#endif
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243#endif
244
245// ************************************************************************* //
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
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
EulerDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
TypeName("Euler")
Runtime type information.
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
EulerDdtScheme(const fvMesh &mesh)
Construct from mesh.
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
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