steadyStateDdtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 "steadyStateDdtScheme.H"
27 #include "fvcDiv.H"
28 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fv
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 tmp<GeometricField<Type, fvPatchField, volMesh>>
45 (
46  const dimensioned<Type>& dt
47 )
48 {
50  (
51  "ddt("+dt.name()+')',
52  mesh(),
54  (
55  "0",
56  dt.dimensions()/dimTime,
57  Zero
58  )
59  );
60 }
61 
62 
63 template<class Type>
66 (
68 )
69 {
71  (
72  "ddt("+vf.name()+')',
73  mesh(),
75  (
76  "0",
77  vf.dimensions()/dimTime,
78  Zero
79  )
80  );
81 }
82 
83 
84 template<class Type>
87 (
88  const dimensionedScalar& rho,
90 )
91 {
93  (
94  "ddt("+rho.name()+','+vf.name()+')',
95  mesh(),
97  (
98  "0",
99  rho.dimensions()*vf.dimensions()/dimTime,
100  Zero
101  )
102  );
103 }
104 
105 
106 template<class Type>
109 (
110  const volScalarField& rho,
112 )
113 {
115  (
116  "ddt("+rho.name()+','+vf.name()+')',
117  mesh(),
119  (
120  "0",
121  rho.dimensions()*vf.dimensions()/dimTime,
122  Zero
123  )
124  );
125 }
126 
127 
128 template<class Type>
131 (
132  const volScalarField& alpha,
133  const volScalarField& rho,
135 )
136 {
138  (
139  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
140  mesh(),
142  (
143  "0",
144  rho.dimensions()*vf.dimensions()/dimTime,
145  Zero
146  )
147  );
148 }
149 
150 
151 template<class Type>
154 (
156 )
157 {
158  return tmp<fvMatrix<Type>>
159  (
160  new fvMatrix<Type>
161  (
162  vf,
163  vf.dimensions()*dimVol/dimTime
164  )
165  );
166 }
167 
168 
169 template<class Type>
172 (
173  const dimensionedScalar& rho,
175 )
176 {
177  return tmp<fvMatrix<Type>>
178  (
179  new fvMatrix<Type>
180  (
181  vf,
182  rho.dimensions()*vf.dimensions()*dimVol/dimTime
183  )
184  );
185 }
186 
187 
188 template<class Type>
191 (
192  const volScalarField& rho,
194 )
195 {
196  return tmp<fvMatrix<Type>>
197  (
198  new fvMatrix<Type>
199  (
200  vf,
201  rho.dimensions()*vf.dimensions()*dimVol/dimTime
202  )
203  );
204 }
205 
206 
207 template<class Type>
210 (
211  const volScalarField& alpha,
212  const volScalarField& rho,
214 )
215 {
216  return tmp<fvMatrix<Type>>
217  (
218  new fvMatrix<Type>
219  (
220  vf,
221  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
222  )
223  );
224 }
225 
226 
227 template<class Type>
230 (
233 )
234 {
235  return fluxFieldType::New
236  (
237  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
238  mesh(),
240  (
241  "0",
243  Zero
244  )
245  );
246 }
247 
248 
249 template<class Type>
252 (
254  const fluxFieldType& phi
255 )
256 {
257  return fluxFieldType::New
258  (
259  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
260  mesh(),
262  (
263  "0",
264  phi.dimensions()/dimTime,
265  Zero
266  )
267  );
268 }
269 
270 
271 template<class Type>
274 (
275  const volScalarField& rho,
278 )
279 {
280  return fluxFieldType::New
281  (
282  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
283  mesh(),
285  (
286  "0",
288  Zero
289  )
290  );
291 }
292 
293 
294 template<class Type>
297 (
298  const volScalarField& rho,
300  const fluxFieldType& phi
301 )
302 {
303  return fluxFieldType::New
304  (
305  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
306  mesh(),
308  (
309  "0",
310  phi.dimensions()/dimTime,
311  Zero
312  )
313  );
314 }
315 
316 
317 template<class Type>
319 (
321 )
322 {
324  (
325  "meshPhi",
326  mesh(),
328  );
329 }
330 
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 } // End namespace fv
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 } // End namespace Foam
339 
340 // ************************************************************************* //
const word & name() const
Return name.
Definition: IOobject.H:303
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
Generic GeometricField class.
Generic dimensioned Type class.
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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:72
static const zero Zero
Definition: zero.H:97
Calculate the divergence of the given field.
const word & name() const
Return const reference to name.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.
const dimensionSet & dimensions() const
Return const reference to dimensions.
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.