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 | Copyright (C) 2011-2016 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<GeometricField<Type, fvPatchField, volMesh>>
44 ddt
45 (
46  const dimensioned<Type> dt,
47  const fvMesh& mesh
48 )
49 {
51  (
52  mesh,
53  mesh.ddtScheme("ddt(" + dt.name() + ')')
54  ).ref().fvcDdt(dt);
55 }
56 
57 
58 template<class Type>
60 ddt
61 (
63 )
64 {
66  (
67  vf.mesh(),
68  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
69  ).ref().fvcDdt(vf);
70 }
71 
72 
73 template<class Type>
75 ddt
76 (
77  const dimensionedScalar& rho,
79 )
80 {
82  (
83  vf.mesh(),
84  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
85  ).ref().fvcDdt(rho, vf);
86 }
87 
88 
89 template<class Type>
91 ddt
92 (
93  const volScalarField& rho,
95 )
96 {
98  (
99  vf.mesh(),
100  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
101  ).ref().fvcDdt(rho, vf);
102 }
103 
104 
105 template<class Type>
107 ddt
108 (
109  const volScalarField& alpha,
110  const volScalarField& rho,
112 )
113 {
115  (
116  vf.mesh(),
117  vf.mesh().ddtScheme
118  (
119  "ddt("
120  + alpha.name() + ','
121  + rho.name() + ','
122  + vf.name() + ')'
123  )
124  ).ref().fvcDdt(alpha, rho, vf);
125 }
126 
127 
128 template<class Type>
130 ddt
131 (
132  const one&,
134 )
135 {
136  return ddt(vf);
137 }
138 
139 
140 template<class Type>
142 ddt
143 (
145  const one&
146 )
147 {
148  return ddt(vf);
149 }
150 
151 
152 template<class Type>
154 ddtCorr
155 (
158 )
159 {
161  (
162  U.mesh(),
163  U.mesh().ddtScheme("ddt(" + U.name() + ')')
164  ).ref().fvcDdtUfCorr(U, Uf);
165 }
166 
167 
168 template<class Type>
170 ddtCorr
171 (
173  const GeometricField
174  <
175  typename flux<Type>::type,
176  fvsPatchField,
178  >& phi
179 )
180 {
182  (
183  U.mesh(),
184  U.mesh().ddtScheme("ddt(" + U.name() + ')')
185  ).ref().fvcDdtPhiCorr(U, phi);
186 }
187 
188 
189 template<class Type>
191 ddtCorr
192 (
193  const volScalarField& rho,
196 )
197 {
199  (
200  U.mesh(),
201  U.mesh().ddtScheme("ddt(" + U.name() + ')')
202  ).ref().fvcDdtUfCorr(rho, U, Uf);
203 }
204 
205 
206 template<class Type>
208 ddtCorr
209 (
210  const volScalarField& rho,
212  const GeometricField
213  <
214  typename flux<Type>::type,
215  fvsPatchField,
217  >& phi
218 )
219 {
221  (
222  U.mesh(),
223  U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
224  ).ref().fvcDdtPhiCorr(rho, U, phi);
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace fvc
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace Foam
235 
236 // ************************************************************************* //
surfaceScalarField & phi
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:47
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
Generic GeometricField class.
Generic dimensioned Type class.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Calculate the first temporal derivative.
const word & name() const
Return const reference to name.
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
const Mesh & mesh() const
Return mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:54
ITstream & ddtScheme(const word &name) const
Definition: fvSchemes.C:360
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank)-2 >::type type
Definition: products.H:115
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
A class representing the concept of 1 (scalar(1.0)) used to avoid unnecessary manipulations for objec...
Definition: one.H:50