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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::fv::boundedDdtScheme
28 
29 Group
30  grpFvDdtSchemes
31 
32 Description
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 
42 SourceFiles
43  boundedDdtScheme.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef boundedDdtScheme_H
48 #define boundedDdtScheme_H
49 
50 #include "ddtScheme.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace fv
60 {
61 
62 /*---------------------------------------------------------------------------*\
63  Class boundedDdtScheme Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 template<class Type>
67 class boundedDdtScheme
68 :
69  public fv::ddtScheme<Type>
70 {
71  // Private data
72 
73  tmp<fv::ddtScheme<Type>> scheme_;
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 
85 public:
86 
87  //- Runtime type information
88  TypeName("bounded");
89 
90 
91  // Constructors
92 
93  //- Construct from mesh and Istream
94  boundedDdtScheme(const fvMesh& mesh, Istream& is)
95  :
96  ddtScheme<Type>(mesh, is),
97  scheme_
98  (
100  )
101  {}
102 
103 
104  // Member Functions
105 
106  //- Return mesh reference
107  const fvMesh& mesh() const
108  {
109  return fv::ddtScheme<Type>::mesh();
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  );
164 
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 
200 template<>
202 (
205 );
206 
207 template<>
209 (
210  const volScalarField& U,
211  const surfaceScalarField& phi
212 );
213 
214 template<>
216 (
217  const volScalarField& rho,
218  const volScalarField& U,
219  const surfaceScalarField& Uf
220 );
221 
222 template<>
224 (
225  const volScalarField& rho,
226  const volScalarField& U,
227  const surfaceScalarField& phi
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 // ************************************************************************* //
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Bounded form of the selected ddt scheme.
const fvMesh & mesh() const
Return mesh reference.
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
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
TypeName("bounded")
Runtime type information.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:81
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:179
labelList fv(nPoints)
ddtScheme< Type >::fluxFieldType fluxFieldType
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
autoPtr< surfaceVectorField > Uf
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const volScalarField & psi
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)