boundedDdtScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fv::boundedDdtScheme
26 
27 Description
28  Bounded form of the selected ddt scheme.
29 
30  Boundedness is achieved by subtracting ddt(phi)*vf or Sp(ddt(rho), vf)
31  which is non-conservative if ddt(rho) != 0 but conservative otherwise.
32 
33  Can be used for the ddt of bounded scalar properties to improve stability
34  if insufficient convergence of the pressure equation causes temporary
35  divergence of the flux field.
36 
37 SourceFiles
38  boundedDdtScheme.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef boundedDdtScheme_H
43 #define boundedDdtScheme_H
44 
45 #include "ddtScheme.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace fv
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class boundedDdtScheme Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class Type>
62 class boundedDdtScheme
63 :
64  public fv::ddtScheme<Type>
65 {
66  // Private Data
67 
68  tmp<fv::ddtScheme<Type>> scheme_;
69 
70 
71 public:
72 
73  //- Runtime type information
74  TypeName("bounded");
75 
76 
77  // Constructors
78 
79  //- Construct from mesh and Istream
80  boundedDdtScheme(const fvMesh& mesh, Istream& is)
81  :
82  ddtScheme<Type>(mesh, is),
83  scheme_
84  (
85  fv::ddtScheme<Type>::New(mesh, is)
86  )
87  {}
88 
89  //- Disallow default bitwise copy construction
90  boundedDdtScheme(const boundedDdtScheme&) = delete;
91 
92 
93  // Member Functions
94 
95  //- Return mesh reference
96  const fvMesh& mesh() const
97  {
99  }
100 
102  (
103  const dimensioned<Type>&
104  );
105 
107  (
109  );
110 
112  (
113  const dimensionedScalar&,
115  );
116 
118  (
119  const volScalarField&,
121  );
122 
124  (
125  const volScalarField& alpha,
126  const volScalarField& rho,
128  );
129 
130  virtual tmp<fvMatrix<Type>> fvmDdt
131  (
133  );
134 
135  virtual tmp<fvMatrix<Type>> fvmDdt
136  (
137  const dimensionedScalar&,
139  );
140 
141  virtual tmp<fvMatrix<Type>> fvmDdt
142  (
143  const volScalarField&,
145  );
146 
147  virtual tmp<fvMatrix<Type>> fvmDdt
148  (
149  const volScalarField& alpha,
150  const volScalarField& rho,
152  );
155 
157  (
160  );
161 
163  (
165  const fluxFieldType& phi
166  );
167 
169  (
170  const volScalarField& rho,
173  );
174 
176  (
177  const volScalarField& rho,
179  const fluxFieldType& phi
180  );
181 
183  (
185  );
186 
187 
188  // Member Operators
189 
190  //- Disallow default bitwise assignment
191  void operator=(const boundedDdtScheme&) = delete;
192 };
193 
194 
195 template<>
197 (
200 );
201 
202 template<>
204 (
205  const volScalarField& U,
206  const surfaceScalarField& phi
207 );
208 
209 template<>
211 (
212  const volScalarField& rho,
213  const volScalarField& U,
214  const surfaceScalarField& Uf
215 );
216 
217 template<>
219 (
220  const volScalarField& rho,
221  const volScalarField& U,
222  const surfaceScalarField& phi
223 );
224 
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 } // End namespace fv
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 } // End namespace Foam
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #ifdef NoRepository
237  #include "boundedDdtScheme.C"
238 #endif
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #endif
243 
244 // ************************************************************************* //
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Bounded form of the selected ddt scheme.
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Generic GeometricField class.
Generic dimensioned Type class.
TypeName("bounded")
Runtime type information.
boundedDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
labelList fv(nPoints)
ddtScheme< Type >::fluxFieldType fluxFieldType
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
phi
Definition: correctPhi.H:3
autoPtr< surfaceVectorField > Uf
void operator=(const boundedDdtScheme &)=delete
Disallow default bitwise assignment.
U
Definition: pEqn.H:72
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)