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 (
133 )
134 {
136  (
137  sf.mesh(),
138  sf.mesh().ddtScheme("ddt(" + sf.name() + ')')
139  ).ref().fvcDdt(sf);
140 }
141 
142 
143 template<class Type>
145 ddt
146 (
147  const one&,
149 )
150 {
151  return ddt(vf);
152 }
153 
154 
155 template<class Type>
157 ddt
158 (
160  const one&
161 )
162 {
163  return ddt(vf);
164 }
165 
166 
167 template<class Type>
169 ddtCorr
170 (
173 )
174 {
176  (
177  U.mesh(),
178  U.mesh().ddtScheme("ddt(" + U.name() + ')')
179  ).ref().fvcDdtUfCorr(U, Uf);
180 }
181 
182 
183 template<class Type>
185 ddtCorr
186 (
188  const GeometricField
189  <
190  typename flux<Type>::type,
191  fvsPatchField,
193  >& phi
194 )
195 {
197  (
198  U.mesh(),
199  U.mesh().ddtScheme("ddt(" + U.name() + ')')
200  ).ref().fvcDdtPhiCorr(U, phi);
201 }
202 
203 
204 template<class Type>
206 ddtCorr
207 (
208  const volScalarField& rho,
211 )
212 {
214  (
215  U.mesh(),
216  U.mesh().ddtScheme("ddt(" + U.name() + ')')
217  ).ref().fvcDdtUfCorr(rho, U, Uf);
218 }
219 
220 
221 template<class Type>
223 ddtCorr
224 (
225  const volScalarField& rho,
227  const GeometricField
228  <
229  typename flux<Type>::type,
230  fvsPatchField,
232  >& phi
233 )
234 {
236  (
237  U.mesh(),
238  U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
239  ).ref().fvcDdtPhiCorr(rho, U, phi);
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace fvc
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 } // End namespace Foam
250 
251 // ************************************************************************* //
surfaceScalarField & phi
const word & name() const
Return name.
Definition: IOobject.H:291
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:170
Generic GeometricField class.
ITstream & ddtScheme(const word &name) const
Definition: fvSchemes.C:360
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 Mesh & mesh() const
Return mesh.
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:115
Namespace for OpenFOAM.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50