backwardDdtScheme.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) 2017 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::backwardDdtScheme
29
30Group
31 grpFvDdtSchemes
32
33Description
34 Second-order backward-differencing ddt using the current and
35 two previous time-step values.
36
37SourceFiles
38 backwardDdtScheme.C
39
40\*---------------------------------------------------------------------------*/
41
42#ifndef backwardDdtScheme_H
43#define backwardDdtScheme_H
44
45#include "ddtScheme.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace fv
55{
56
57/*---------------------------------------------------------------------------*\
58 Class backwardDdtScheme Declaration
59\*---------------------------------------------------------------------------*/
60
61template<class Type>
63:
64 public fv::ddtScheme<Type>
65{
66 // Private Member Functions
67
68 //- Return the current time-step
69 scalar deltaT_() const;
70
71 //- Return the previous time-step
72 scalar deltaT0_() const;
73
74 //- Return the previous time-step or GREAT if the old timestep field
75 // wasn't available in which case Euler ddt is used
76 template<class GeoField>
77 scalar deltaT0_(const GeoField&) const;
78
79 //- No copy construct
80 backwardDdtScheme(const backwardDdtScheme&) = delete;
81
82 //- No copy assignment
83 void operator=(const backwardDdtScheme&) = delete;
84
85
86public:
87
88 //- Runtime type information
89 TypeName("backward");
90
91
92 // Constructors
93
94 //- Construct from mesh
96 :
97 ddtScheme<Type>(mesh)
98 {}
99
100 //- Construct from mesh and Istream
102 :
103 ddtScheme<Type>(mesh, is)
104 {
105 if (is.good() && !is.eof())
106 {
107 is >> this->ddtPhiCoeff_;
108 }
109
110 // Ensure the old-old-time cell volumes are available
111 // for moving meshes
112 if (mesh.moving())
113 {
114 mesh.V00();
115 }
116 }
117
118
119 // Member Functions
120
121 //- Return mesh reference
122 const fvMesh& mesh() const
123 {
125 }
126
128 (
129 const dimensioned<Type>&
130 );
131
133 (
135 );
136
138 (
139 const dimensionedScalar&,
141 );
142
144 (
145 const volScalarField&,
147 );
148
150 (
151 const volScalarField& alpha,
152 const volScalarField& rho,
154 );
155
157 (
159 );
160
162 (
163 const dimensionedScalar&,
165 );
166
168 (
169 const volScalarField&,
171 );
172
174 (
175 const volScalarField& alpha,
176 const volScalarField& rho,
178 );
181
183 (
186 );
187
189 (
191 const fluxFieldType& phi
192 );
193
195 (
196 const volScalarField& rho,
199 );
200
202 (
203 const volScalarField& rho,
205 const fluxFieldType& phi
206 );
207
209 (
211 );
212};
213
214
215template<>
217(
220);
221
222template<>
224(
225 const volScalarField& U,
227);
228
229template<>
231(
232 const volScalarField& rho,
233 const volScalarField& U,
235);
236
237template<>
239(
240 const volScalarField& rho,
241 const volScalarField& U,
243);
244
245
246// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247
248} // End namespace fv
249
250// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251
252} // End namespace Foam
253
254// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255
256#ifdef NoRepository
257 #include "backwardDdtScheme.C"
258#endif
259
260// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261
262#endif
263
264// ************************************************************************* //
surfaceScalarField & phi
Generic GeometricField class.
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
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
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Second-order backward-differencing ddt using the current and two previous time-step values.
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("backward")
Runtime type information.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
backwardDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
const fvMesh & mesh() const
Return mesh reference.
backwardDdtScheme(const fvMesh &mesh)
Construct from mesh.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Abstract base class for ddt schemes.
Definition: ddtScheme.H:87
scalar ddtPhiCoeff_
Input for fvcDdtPhiCoeff.
Definition: ddtScheme.H:98
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:534
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