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-2023 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 
101  virtual tmp<VolField<Type>> fvcDdt
102  (
103  const dimensioned<Type>&
104  );
105 
106  virtual tmp<VolField<Type>> fvcDdt
107  (
108  const VolField<Type>&
109  );
110 
111  virtual tmp<VolField<Type>> fvcDdt
112  (
113  const dimensionedScalar&,
114  const VolField<Type>&
115  );
116 
117  virtual tmp<VolField<Type>> fvcDdt
118  (
119  const volScalarField&,
120  const VolField<Type>&
121  );
122 
123  virtual tmp<VolField<Type>> fvcDdt
124  (
125  const volScalarField& alpha,
126  const volScalarField& rho,
127  const VolField<Type>& psi
128  );
129 
130  virtual tmp<fvMatrix<Type>> fvmDdt
131  (
132  const VolField<Type>&
133  );
134 
135  virtual tmp<fvMatrix<Type>> fvmDdt
136  (
137  const dimensionedScalar&,
138  const VolField<Type>&
139  );
140 
141  virtual tmp<fvMatrix<Type>> fvmDdt
142  (
143  const volScalarField&,
144  const VolField<Type>&
145  );
146 
147  virtual tmp<fvMatrix<Type>> fvmDdt
148  (
149  const volScalarField& alpha,
150  const volScalarField& rho,
151  const VolField<Type>& psi
152  );
153 
155 
157  (
158  const VolField<Type>& U,
159  const SurfaceField<Type>& Uf
160  );
161 
163  (
164  const VolField<Type>& U,
165  const fluxFieldType& phi
166  );
167 
169  (
170  const volScalarField& rho,
171  const VolField<Type>& U,
172  const SurfaceField<Type>& rhoUf
173  );
174 
176  (
177  const volScalarField& rho,
178  const VolField<Type>& U,
179  const fluxFieldType& phi
180  );
181 
183  (
184  const volScalarField& alpha,
185  const volScalarField& rho,
186  const VolField<Type>& U,
187  const SurfaceField<Type>& Uf
188  )
189  {
191  return fluxFieldType::null();
192  }
193 
195  (
196  const volScalarField& alpha,
197  const volScalarField& rho,
198  const VolField<Type>& U,
199  const fluxFieldType& phi
200  )
201  {
203  return fluxFieldType::null();
204  }
205 
207  (
208  const VolField<Type>&
209  );
210 
211  virtual tmp<scalarField> meshPhi
212  (
213  const VolField<Type>&,
214  const label patchi
215  );
216 
217 
218  // Member Operators
219 
220  //- Disallow default bitwise assignment
221  void operator=(const boundedDdtScheme&) = delete;
222 };
223 
224 
225 template<>
227 (
228  const VolField<scalar>& U,
229  const SurfaceField<scalar>& Uf
230 );
231 
232 template<>
234 (
235  const volScalarField& U,
236  const surfaceScalarField& phi
237 );
238 
239 template<>
241 (
242  const volScalarField& rho,
243  const volScalarField& U,
244  const surfaceScalarField& Uf
245 );
246 
247 template<>
249 (
250  const volScalarField& rho,
251  const volScalarField& U,
252  const surfaceScalarField& phi
253 );
254 
255 template<>
257 (
258  const volScalarField& alpha,
259  const volScalarField& rho,
260  const volScalarField& U,
261  const surfaceScalarField& Uf
262 );
263 
264 template<>
266 (
267  const volScalarField& alpha,
268  const volScalarField& rho,
269  const volScalarField& U,
270  const surfaceScalarField& phi
271 );
272 
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 } // End namespace fv
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 } // End namespace Foam
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 #ifdef NoRepository
285  #include "boundedDdtScheme.C"
286 #endif
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 #endif
291 
292 // ************************************************************************* //
Generic GeometricField class.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Bounded form of the selected ddt scheme.
boundedDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
void operator=(const boundedDdtScheme &)=delete
Disallow default bitwise assignment.
TypeName("bounded")
Runtime type information.
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
Abstract base class for ddt schemes.
Definition: ddtScheme.H:67
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:45
A class for managing temporary objects.
Definition: tmp.H:55
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
label patchi
const volScalarField & psi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
labelList fv(nPoints)