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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "steadyStateDdtScheme.H"
30 #include "fvcDiv.H"
31 #include "fvMatrices.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace fv
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<class Type>
48 (
49  const dimensioned<Type>& dt
50 )
51 {
53  (
54  IOobject
55  (
56  "ddt("+dt.name()+')',
57  mesh().time().timeName(),
58  mesh().thisDb()
59  ),
60  mesh(),
62  );
63 }
64 
65 
66 template<class Type>
69 (
71 )
72 {
74  (
75  IOobject
76  (
77  "ddt("+vf.name()+')',
78  mesh().time().timeName(),
79  mesh().thisDb()
80  ),
81  mesh(),
83  );
84 }
85 
86 
87 template<class Type>
90 (
91  const dimensionedScalar& rho,
93 )
94 {
96  (
97  IOobject
98  (
99  "ddt("+rho.name()+','+vf.name()+')',
100  mesh().time().timeName(),
101  mesh().thisDb()
102  ),
103  mesh(),
104  dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
105  );
106 }
107 
108 
109 template<class Type>
112 (
113  const volScalarField& rho,
115 )
116 {
118  (
119  IOobject
120  (
121  "ddt("+rho.name()+','+vf.name()+')',
122  mesh().time().timeName(),
123  mesh().thisDb()
124  ),
125  mesh(),
126  dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
127  );
128 }
129 
130 
131 template<class Type>
134 (
135  const volScalarField& alpha,
136  const volScalarField& rho,
138 )
139 {
141  (
142  IOobject
143  (
144  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
145  mesh().time().timeName(),
146  mesh().thisDb()
147  ),
148  mesh(),
149  dimensioned<Type>(rho.dimensions()*vf.dimensions()/dimTime, Zero)
150  );
151 }
152 
153 
154 template<class Type>
157 (
159 )
160 {
161  tmp<fvMatrix<Type>> tfvm
162  (
163  new fvMatrix<Type>
164  (
165  vf,
167  )
168  );
169 
170  return tfvm;
171 }
172 
173 
174 template<class Type>
177 (
178  const dimensionedScalar& rho,
180 )
181 {
182  tmp<fvMatrix<Type>> tfvm
183  (
184  new fvMatrix<Type>
185  (
186  vf,
187  rho.dimensions()*vf.dimensions()*dimVol/dimTime
188  )
189  );
190 
191  return tfvm;
192 }
193 
194 
195 template<class Type>
198 (
199  const volScalarField& rho,
201 )
202 {
203  tmp<fvMatrix<Type>> tfvm
204  (
205  new fvMatrix<Type>
206  (
207  vf,
208  rho.dimensions()*vf.dimensions()*dimVol/dimTime
209  )
210  );
211 
212  return tfvm;
213 }
214 
215 
216 template<class Type>
219 (
220  const volScalarField& alpha,
221  const volScalarField& rho,
223 )
224 {
225  tmp<fvMatrix<Type>> tfvm
226  (
227  new fvMatrix<Type>
228  (
229  vf,
230  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
231  )
232  );
233 
234  return tfvm;
235 }
236 
237 
238 template<class Type>
241 (
244 )
245 {
246  tmp<fluxFieldType> tCorr
247  (
248  new fluxFieldType
249  (
250  IOobject
251  (
252  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
253  mesh().time().timeName(),
254  mesh().thisDb()
255  ),
256  mesh(),
257  Foam::zero{}, // value
258  Uf.dimensions()*dimArea/dimTime
259  )
260  );
261 
262  tCorr.ref().setOriented();
263 
264  return tCorr;
265 }
266 
267 
268 template<class Type>
271 (
273  const fluxFieldType& phi
274 )
275 {
276  tmp<fluxFieldType> tCorr
277  (
278  new fluxFieldType
279  (
280  IOobject
281  (
282  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
283  mesh().time().timeName(),
284  mesh().thisDb()
285  ),
286  mesh(),
287  Foam::zero{}, // value
288  phi.dimensions()/dimTime
289  )
290  );
291 
292  tCorr.ref().setOriented();
293 
294  return tCorr;
295 }
296 
297 
298 template<class Type>
301 (
302  const volScalarField& rho,
305 )
306 {
307  tmp<fluxFieldType> tCorr
308  (
309  new fluxFieldType
310  (
311  IOobject
312  (
313  "ddtCorr("
314  + rho.name()
315  + ',' + U.name() + ',' + Uf.name() + ')',
316  mesh().time().timeName(),
317  mesh().thisDb()
318  ),
319  mesh(),
320  Foam::zero{}, // value
321  Uf.dimensions()*dimArea/dimTime
322  )
323  );
324 
325  tCorr.ref().setOriented();
326 
327  return tCorr;
328 }
329 
330 
331 template<class Type>
334 (
335  const volScalarField& rho,
337  const fluxFieldType& phi
338 )
339 {
340  tmp<fluxFieldType> tCorr
341  (
342  new fluxFieldType
343  (
344  IOobject
345  (
346  "ddtCorr("
347  + rho.name()
348  + ',' + U.name() + ',' + phi.name() + ')',
349  mesh().time().timeName(),
350  mesh().thisDb()
351  ),
352  mesh(),
353  Foam::zero{}, // value
354  phi.dimensions()/dimTime
355  )
356  );
357 
358  tCorr.ref().setOriented();
360  return tCorr;
361 }
362 
363 
364 template<class Type>
366 (
368 )
369 {
370  auto tphi
371  (
373  (
374  IOobject
375  (
376  "meshPhi",
377  mesh().time().timeName(),
378  mesh().thisDb(),
382  ),
383  mesh(),
385  )
386  );
387  tphi.ref().setOriented();
388  return tphi;
389 }
390 
391 
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393 
394 } // End namespace fv
395 
396 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
397 
398 } // End namespace Foam
399 
400 // ************************************************************************* //
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
Generic GeometricField class.
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
labelList fv(nPoints)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
autoPtr< surfaceVectorField > Uf
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
Calculate the divergence of the given field.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
const word & name() const noexcept
Return const reference to name.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127