fvcDdt.C
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) 2011-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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvcDdt.H"
27 #include "fvMesh.H"
28 #include "ddtScheme.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fvc
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 tmp<VolField<Type>>
45 (
46  const dimensioned<Type> dt,
47  const fvMesh& mesh
48 )
49 {
51  (
52  mesh,
53  mesh.schemes().ddt("ddt(" + dt.name() + ')')
54  ).ref().fvcDdt(dt);
55 }
56 
57 
58 template<class Type>
61 (
62  const VolField<Type>& vf
63 )
64 {
66  (
67  vf.mesh(),
68  vf.mesh().schemes().ddt("ddt(" + vf.name() + ')')
69  ).ref().fvcDdt(vf);
70 }
71 
72 
73 template<class Type>
76 (
77  const dimensionedScalar& rho,
78  const VolField<Type>& vf
79 )
80 {
82  (
83  vf.mesh(),
84  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
85  ).ref().fvcDdt(rho, vf);
86 }
87 
88 
89 template<class Type>
92 (
93  const volScalarField& rho,
94  const VolField<Type>& vf
95 )
96 {
98  (
99  vf.mesh(),
100  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
101  ).ref().fvcDdt(rho, vf);
102 }
103 
104 
105 template<class Type>
108 (
109  const one&,
110  const VolField<Type>& vf
111 )
112 {
113  return ddt(vf);
114 }
115 
116 
117 template<class Type>
120 (
121  const volScalarField& alpha,
122  const volScalarField& rho,
123  const VolField<Type>& vf
124 )
125 {
127  (
128  vf.mesh(),
129  vf.mesh().schemes().ddt
130  (
131  "ddt("
132  + alpha.name() + ','
133  + rho.name() + ','
134  + vf.name() + ')'
135  )
136  ).ref().fvcDdt(alpha, rho, vf);
137 }
138 
139 
140 template<class Type>
143 (
144  const one&,
145  const one&,
146  const VolField<Type>& vf
147 )
148 {
149  return ddt(vf);
150 }
151 
152 
153 template<class Type>
156 (
157  const one&,
158  const volScalarField& rho,
159  const VolField<Type>& vf
160 )
161 {
162  return ddt(rho, vf);
163 }
164 
165 
166 template<class Type>
169 (
170  const volScalarField& alpha,
171  const one&,
172  const VolField<Type>& vf
173 )
174 {
175  return ddt(alpha, vf);
176 }
177 
178 
179 template<class Type>
182 (
183  const SurfaceField<Type>& sf
184 )
185 {
187  (
188  sf.mesh(),
189  sf.mesh().schemes().ddt("ddt(" + sf.name() + ')')
190  ).ref().fvcDdt(sf);
191 }
192 
193 
194 template<class Type>
196 (
197  const VolField<Type>& U,
198  const SurfaceField<Type>& Uf
199 )
200 {
202  (
203  U.mesh(),
204  U.mesh().schemes().ddt("ddt(" + U.name() + ')')
205  ).ref().fvcDdtUfCorr(U, Uf);
206 }
207 
208 
209 template<class Type>
211 (
212  const VolField<Type>& U,
213  const SurfaceField<typename Foam::flux<Type>::type>& phi
214 )
215 {
217  (
218  U.mesh(),
219  U.mesh().schemes().ddt("ddt(" + U.name() + ')')
220  ).ref().fvcDdtPhiCorr(U, phi);
221 }
222 
223 
224 template<class Type>
226 (
227  const VolField<Type>& U,
228  const SurfaceField<typename Foam::flux<Type>::type>& phi,
229  const autoPtr<SurfaceField<Type>>& Uf
230 )
231 {
232  if (Uf.valid())
233  {
234  return ddtCorr(U, Uf());
235  }
236  else
237  {
238  return ddtCorr(U, phi);
239  }
240 }
241 
242 
243 template<class Type>
245 (
246  const volScalarField& rho,
247  const VolField<Type>& U,
248  const SurfaceField<Type>& rhoUf
249 )
250 {
252  (
253  U.mesh(),
254  U.mesh().schemes().ddt
255  (
256  "ddt(" + rho.name() + U.name() + ')'
257  )
258  ).ref().fvcDdtUfCorr(rho, U, rhoUf);
259 }
260 
261 
262 template<class Type>
264 (
265  const volScalarField& rho,
266  const VolField<Type>& U,
267  const SurfaceField<typename Foam::flux<Type>::type>& phi
268 )
269 {
271  (
272  U.mesh(),
273  U.mesh().schemes().ddt("ddt(" + rho.name() + ',' + U.name() + ')')
274  ).ref().fvcDdtPhiCorr(rho, U, phi);
275 }
276 
277 
278 template<class Type>
280 (
281  const volScalarField& rho,
282  const VolField<Type>& U,
283  const SurfaceField<typename Foam::flux<Type>::type>& phi,
284  const autoPtr<SurfaceField<Type>>& rhoUf
285 )
286 {
287  if (rhoUf.valid())
288  {
289  return ddtCorr(rho, U, rhoUf());
290  }
291  else
292  {
293  return ddtCorr(rho, U, phi);
294  }
295 }
296 
297 
298 template<class Type>
300 (
301  const volScalarField& alpha,
302  const volScalarField& rho,
303  const VolField<Type>& U,
304  const SurfaceField<Type>& Uf
305 )
306 {
308  (
309  U.mesh(),
310  U.mesh().schemes().ddt
311  (
312  "ddt(" + alpha.name() + rho.name() + ',' + U.name() + ')'
313  )
314  ).ref().fvcDdtUfCorr(alpha, rho, U, Uf);
315 }
316 
317 
318 template<class Type>
320 (
321  const volScalarField& alpha,
322  const volScalarField& rho,
323  const VolField<Type>& U,
324  const SurfaceField<typename Foam::flux<Type>::type>& phi
325 )
326 {
328  (
329  U.mesh(),
330  U.mesh().schemes().ddt
331  (
332  "ddt(" + alpha.name() + rho.name() + ',' + U.name() + ')'
333  )
334  ).ref().fvcDdtPhiCorr(alpha, rho, U, phi);
335 }
336 
337 
338 template<class Type>
340 (
341  const volScalarField& alpha,
342  const volScalarField& rho,
343  const VolField<Type>& U,
344  const SurfaceField<typename Foam::flux<Type>::type>& phi,
345  const autoPtr<SurfaceField<Type>>& Uf
346 )
347 {
348  if (Uf.valid())
349  {
350  return ddtCorr(alpha, rho, U, Uf());
351  }
352  else
353  {
354  return ddtCorr(alpha, rho, U, phi);
355  }
356 }
357 
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 } // End namespace fvc
362 
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 
365 } // End namespace Foam
366 
367 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Generic dimensioned Type class.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:456
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:45
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:115
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:51
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the first temporal derivative.
volScalarField sf(fieldObject, mesh)
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
Namespace for OpenFOAM.