boundedDdtScheme.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) 2012-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::boundedDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 Bounded form of the selected ddt scheme.
34
35 Boundedness is achieved by subtracting ddt(phi)*vf or Sp(ddt(rho), vf)
36 which is non-conservative if ddt(rho) != 0 but conservative otherwise.
37
38 Can be used for the ddt of bounded scalar properties to improve stability
39 if insufficient convergence of the pressure equation causes temporary
40 divergence of the flux field.
41
42SourceFiles
43 boundedDdtScheme.C
44
45\*---------------------------------------------------------------------------*/
46
47#ifndef boundedDdtScheme_H
48#define boundedDdtScheme_H
49
50#include "ddtScheme.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace fv
60{
61
62/*---------------------------------------------------------------------------*\
63 Class boundedDdtScheme Declaration
64\*---------------------------------------------------------------------------*/
65
66template<class Type>
68:
69 public fv::ddtScheme<Type>
70{
71 // Private data
72
74
75
76 // Private Member Functions
77
78 //- No copy construct
79 boundedDdtScheme(const boundedDdtScheme&) = delete;
80
81 //- No copy assignment
82 void operator=(const boundedDdtScheme&) = delete;
83
84
85public:
86
87 //- Runtime type information
88 TypeName("bounded");
89
90
91 // Constructors
92
93 //- Construct from mesh and Istream
95 :
96 ddtScheme<Type>(mesh, is),
97 scheme_
98 (
99 fv::ddtScheme<Type>::New(mesh, is)
100 )
101 {}
102
103
104 // Member Functions
105
106 //- Return mesh reference
107 const fvMesh& mesh() const
108 {
110 }
111
113 (
114 const dimensioned<Type>&
115 );
116
118 (
120 );
121
123 (
124 const dimensionedScalar&,
126 );
127
129 (
130 const volScalarField&,
132 );
133
135 (
136 const volScalarField& alpha,
137 const volScalarField& rho,
139 );
140
142 (
144 );
145
147 (
148 const dimensionedScalar&,
150 );
151
153 (
154 const volScalarField&,
156 );
157
159 (
160 const volScalarField& alpha,
161 const volScalarField& rho,
163 );
166
168 (
171 );
172
174 (
176 const fluxFieldType& phi
177 );
178
180 (
181 const volScalarField& rho,
184 );
185
187 (
188 const volScalarField& rho,
190 const fluxFieldType& phi
191 );
192
194 (
196 );
197};
198
199
200template<>
202(
205);
206
207template<>
209(
210 const volScalarField& U,
212);
213
214template<>
216(
217 const volScalarField& rho,
218 const volScalarField& U,
220);
221
222template<>
224(
225 const volScalarField& rho,
226 const volScalarField& U,
228);
229
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace fv
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237} // End namespace Foam
238
239// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240
241#ifdef NoRepository
242 #include "boundedDdtScheme.C"
243#endif
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247#endif
248
249// ************************************************************************* //
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
Bounded form of the selected ddt scheme.
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
boundedDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
TypeName("bounded")
Runtime type information.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
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
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:50
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