steadyStateDdtScheme.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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
29#include "fvcDiv.H"
30#include "fvMatrices.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace fv
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44template<class Type>
45tmp<GeometricField<Type, fvPatchField, volMesh>>
47(
48 const dimensioned<Type>& dt
49)
50{
52 (
54 (
55 "ddt("+dt.name()+')',
56 mesh().time().timeName(),
57 mesh()
58 ),
59 mesh(),
61 );
62}
63
64
65template<class Type>
68(
70)
71{
73 (
75 (
76 "ddt("+vf.name()+')',
77 mesh().time().timeName(),
78 mesh()
79 ),
80 mesh(),
82 );
83}
84
85
86template<class Type>
89(
92)
93{
95 (
97 (
98 "ddt("+rho.name()+','+vf.name()+')',
99 mesh().time().timeName(),
100 mesh()
101 ),
102 mesh(),
103 dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
104 );
105}
106
107
108template<class Type>
111(
112 const volScalarField& rho,
114)
115{
117 (
119 (
120 "ddt("+rho.name()+','+vf.name()+')',
121 mesh().time().timeName(),
122 mesh()
123 ),
124 mesh(),
125 dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
126 );
127}
128
129
130template<class Type>
133(
134 const volScalarField& alpha,
135 const volScalarField& rho,
137)
138{
140 (
142 (
143 "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
144 mesh().time().timeName(),
145 mesh()
146 ),
147 mesh(),
148 dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
149 );
150}
151
152
153template<class Type>
156(
158)
159{
161 (
163 (
164 vf,
166 )
167 );
168
169 return tfvm;
170}
171
172
173template<class Type>
176(
177 const dimensionedScalar& rho,
179)
180{
182 (
184 (
185 vf,
186 rho.dimensions()*vf.dimensions()*dimVol/dimTime
187 )
188 );
189
190 return tfvm;
191}
192
193
194template<class Type>
197(
198 const volScalarField& rho,
200)
201{
203 (
205 (
206 vf,
207 rho.dimensions()*vf.dimensions()*dimVol/dimTime
208 )
209 );
210
211 return tfvm;
212}
213
214
215template<class Type>
218(
219 const volScalarField& alpha,
220 const volScalarField& rho,
222)
223{
225 (
227 (
228 vf,
229 alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
230 )
231 );
232
233 return tfvm;
234}
235
236
237template<class Type>
240(
243)
244{
246 (
247 new fluxFieldType
248 (
250 (
251 "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
252 mesh().time().timeName(),
253 mesh()
254 ),
255 mesh(),
257 (
258 Uf.dimensions()*dimArea/dimTime, Zero
259 )
260 )
261 );
262
263 tCorr.ref().setOriented();
264
265 return tCorr;
266}
267
268
269template<class Type>
272(
274 const fluxFieldType& phi
275)
276{
278 (
279 new fluxFieldType
280 (
282 (
283 "ddtCorr(" + U.name() + ',' + phi.name() + ')',
284 mesh().time().timeName(),
285 mesh()
286 ),
287 mesh(),
289 (
290 phi.dimensions()/dimTime, Zero
291 )
292 )
293 );
294
295 tCorr.ref().setOriented();
296
297 return tCorr;
298}
299
300
301template<class Type>
304(
305 const volScalarField& rho,
308)
309{
311 (
312 new fluxFieldType
313 (
315 (
316 "ddtCorr("
317 + rho.name()
318 + ',' + U.name() + ',' + Uf.name() + ')',
319 mesh().time().timeName(),
320 mesh()
321 ),
322 mesh(),
324 (
325 Uf.dimensions()*dimArea/dimTime, Zero
326 )
327 )
328 );
329
330 tCorr.ref().setOriented();
331
332 return tCorr;
333}
334
335
336template<class Type>
339(
340 const volScalarField& rho,
342 const fluxFieldType& phi
343)
344{
346 (
347 new fluxFieldType
348 (
350 (
351 "ddtCorr("
352 + rho.name()
353 + ',' + U.name() + ',' + phi.name() + ')',
354 mesh().time().timeName(),
355 mesh()
356 ),
357 mesh(),
359 (
360 phi.dimensions()/dimTime, Zero
361 )
362 )
363 );
364
365 tCorr.ref().setOriented();
366
367 return tCorr;
368}
369
370
371template<class Type>
373(
375)
376{
378 (
380 (
381 "meshPhi",
382 mesh().time().timeName(),
383 mesh(),
386 false
387 ),
388 mesh(),
390 );
391}
392
393
394// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395
396} // End namespace fv
397
398// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399
400} // End namespace Foam
401
402// ************************************************************************* //
surfaceScalarField & phi
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
type
Volume classification types.
Definition: volumeType.H:66
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
labelList fv(nPoints)
volScalarField & alpha