steadyStateDdtScheme.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::steadyStateDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 SteadyState implicit/explicit ddt which returns 0.
34
35SourceFiles
36 steadyStateDdtScheme.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef steadyStateDdtScheme_H
41#define steadyStateDdtScheme_H
42
43#include "ddtScheme.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52namespace fv
53{
54
55/*---------------------------------------------------------------------------*\
56 Class steadyStateDdtScheme Declaration
57\*---------------------------------------------------------------------------*/
58
59template<class Type>
61:
62 public fv::ddtScheme<Type>
63{
64 // Private Member Functions
65
66 //- No copy construct
68
69 //- No copy assignment
70 void operator=(const steadyStateDdtScheme&) = delete;
71
72
73public:
74
75 //- Runtime type information
76 TypeName("steadyState");
77
78
79 // Constructors
80
81 //- Construct from mesh
83 :
84 ddtScheme<Type>(mesh)
85 {}
86
87 //- Construct from mesh and Istream
89 :
90 ddtScheme<Type>(mesh, is)
91 {}
92
93
94 // Member Functions
95
96 //- Return mesh reference
97 const fvMesh& mesh() const
98 {
100 }
101
103 (
104 const dimensioned<Type>&
105 );
106
108 (
110 );
111
113 (
114 const dimensionedScalar&,
116 );
117
119 (
120 const volScalarField&,
122 );
123
125 (
126 const volScalarField& alpha,
127 const volScalarField& rho,
129 );
130
132 (
134 );
135
137 (
138 const dimensionedScalar&,
140 );
141
143 (
144 const volScalarField&,
146 );
147
149 (
150 const volScalarField& alpha,
151 const volScalarField& rho,
153 );
156
158 (
161 );
162
164 (
166 const fluxFieldType& phi
167 );
168
170 (
171 const volScalarField& rho,
174 );
175
177 (
178 const volScalarField& rho,
180 const fluxFieldType& phi
181 );
182
184 (
186 );
187};
188
189
190template<>
192(
195);
196
197template<>
199(
200 const volScalarField& U,
202);
203
204template<>
206(
207 const volScalarField& rho,
208 const volScalarField& U,
210);
211
212template<>
214(
215 const volScalarField& rho,
216 const volScalarField& U,
218);
219
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223} // End namespace fv
224
225// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226
227} // End namespace Foam
228
229// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230
231#ifdef NoRepository
232 #include "steadyStateDdtScheme.C"
233#endif
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237#endif
238
239// ************************************************************************* //
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
Abstract base class for ddt schemes.
Definition: ddtScheme.H:87
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
SteadyState implicit/explicit ddt which returns 0.
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
steadyStateDdtScheme(const fvMesh &mesh)
Construct from mesh.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
TypeName("steadyState")
Runtime type information.
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
steadyStateDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
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