fvcDdt.C
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) 2011-2017 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 \*---------------------------------------------------------------------------*/
27 
28 #include "fvcDdt.H"
29 #include "fvMesh.H"
30 #include "ddtScheme.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fvc
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 tmp<GeometricField<Type, fvPatchField, volMesh>>
46 ddt
47 (
48  const dimensioned<Type> dt,
49  const fvMesh& mesh
50 )
51 {
53  (
54  mesh,
55  mesh.ddtScheme("ddt(" + dt.name() + ')')
56  ).ref().fvcDdt(dt);
57 }
58 
59 
60 template<class Type>
62 ddt
63 (
65 )
66 {
68  (
69  vf.mesh(),
70  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
71  ).ref().fvcDdt(vf);
72 }
73 
74 
75 template<class Type>
77 ddt
78 (
79  const dimensionedScalar& rho,
81 )
82 {
84  (
85  vf.mesh(),
86  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
87  ).ref().fvcDdt(rho, vf);
88 }
89 
90 
91 template<class Type>
93 ddt
94 (
95  const volScalarField& rho,
97 )
98 {
100  (
101  vf.mesh(),
102  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
103  ).ref().fvcDdt(rho, vf);
104 }
105 
106 
107 template<class Type>
109 ddt
110 (
111  const volScalarField& alpha,
112  const volScalarField& rho,
114 )
115 {
117  (
118  vf.mesh(),
119  vf.mesh().ddtScheme
120  (
121  "ddt("
122  + alpha.name() + ','
123  + rho.name() + ','
124  + vf.name() + ')'
125  )
126  ).ref().fvcDdt(alpha, rho, vf);
127 }
128 
129 
130 template<class Type>
132 ddt
133 (
135 )
136 {
138  (
139  sf.mesh(),
140  sf.mesh().ddtScheme("ddt(" + sf.name() + ')')
141  ).ref().fvcDdt(sf);
142 }
143 
144 
145 template<class Type>
147 ddt
148 (
149  const one&,
151 )
152 {
153  return ddt(vf);
154 }
155 
156 
157 template<class Type>
159 ddt
160 (
162  const one&
163 )
164 {
165  return ddt(vf);
166 }
167 
168 
169 template<class Type>
171 ddtCorr
172 (
175 )
176 {
178  (
179  U.mesh(),
180  U.mesh().ddtScheme("ddt(" + U.name() + ')')
181  ).ref().fvcDdtUfCorr(U, Uf);
182 }
183 
184 
185 template<class Type>
187 ddtCorr
188 (
190  const GeometricField
191  <
192  typename flux<Type>::type,
195  >& phi
196 )
197 {
199  (
200  U.mesh(),
201  U.mesh().ddtScheme("ddt(" + U.name() + ')')
202  ).ref().fvcDdtPhiCorr(U, phi);
203 }
204 
205 
206 template<class Type>
208 ddtCorr
209 (
211  const GeometricField
212  <
213  typename flux<Type>::type,
216  >& phi,
218 )
219 {
220  if (U.mesh().dynamic())
221  {
222  return ddtCorr(U, Uf());
223  }
224  else
225  {
226  return ddtCorr(U, phi);
227  }
228 }
229 
230 
231 template<class Type>
233 ddtCorr
234 (
235  const volScalarField& rho,
238 )
239 {
241  (
242  U.mesh(),
243  U.mesh().ddtScheme("ddt(" + U.name() + ')')
244  ).ref().fvcDdtUfCorr(rho, U, Uf);
245 }
246 
247 
248 template<class Type>
250 ddtCorr
251 (
252  const volScalarField& rho,
254  const GeometricField
255  <
256  typename flux<Type>::type,
259  >& phi
260 )
261 {
263  (
264  U.mesh(),
265  U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
266  ).ref().fvcDdtPhiCorr(rho, U, phi);
267 }
268 
269 
270 template<class Type>
272 ddtCorr
273 (
274  const volScalarField& rho,
276  const GeometricField
277  <
278  typename flux<Type>::type,
281  >& phi,
283 )
284 {
285  if (U.mesh().dynamic())
286  {
287  return ddtCorr(rho, U, Uf());
288  }
289  else
290  {
291  return ddtCorr(rho, U, phi);
292  }
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace fvc
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace Foam
303 
304 // ************************************************************************* //
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::surfaceMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:49
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
rho
rho
Definition: readInitialConditions.H:88
Foam::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:172
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
ddtScheme.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
fvcDdt.H
Calculate the first temporal derivative.
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::ddtScheme::New
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:50
Foam::innerProduct< vector, Type >::type
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:149