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-2023 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<VolField<Type>>
45 (
46  const dimensioned<Type>& dt
47 )
48 {
49  return VolField<Type>::New
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 (
67  const VolField<Type>& vf
68 )
69 {
70  return VolField<Type>::New
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,
89  const VolField<Type>& vf
90 )
91 {
92  return VolField<Type>::New
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,
111  const VolField<Type>& vf
112 )
113 {
114  return VolField<Type>::New
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,
134  const VolField<Type>& vf
135 )
136 {
137  return VolField<Type>::New
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 (
155  const VolField<Type>& vf
156 )
157 {
158  return tmp<fvMatrix<Type>>
159  (
160  new fvMatrix<Type>
161  (
162  vf,
164  )
165  );
166 }
167 
168 
169 template<class Type>
172 (
173  const dimensionedScalar& rho,
174  const VolField<Type>& vf
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,
193  const VolField<Type>& vf
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,
213  const VolField<Type>& vf
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 (
231  const VolField<Type>& U,
232  const SurfaceField<Type>& Uf
233 )
234 {
235  return fluxFieldType::New
236  (
237  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
238  mesh(),
239  dimensioned<typename flux<Type>::type>
240  (
241  "0",
243  Zero
244  )
245  );
246 }
247 
248 
249 template<class Type>
252 (
253  const VolField<Type>& U,
254  const fluxFieldType& phi
255 )
256 {
257  return fluxFieldType::New
258  (
259  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
260  mesh(),
261  dimensioned<typename flux<Type>::type>
262  (
263  "0",
264  phi.dimensions()/dimTime,
265  Zero
266  )
267  );
268 }
269 
270 
271 template<class Type>
274 (
275  const volScalarField& rho,
276  const VolField<Type>& U,
277  const SurfaceField<Type>& rhoUf
278 )
279 {
280  return fluxFieldType::New
281  (
282  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + rhoUf.name() + ')',
283  mesh(),
284  dimensioned<typename flux<Type>::type>
285  (
286  "0",
287  rhoUf.dimensions()*dimArea/dimTime,
288  Zero
289  )
290  );
291 }
292 
293 
294 template<class Type>
297 (
298  const volScalarField& rho,
299  const VolField<Type>& U,
300  const fluxFieldType& phi
301 )
302 {
303  return fluxFieldType::New
304  (
305  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
306  mesh(),
307  dimensioned<typename flux<Type>::type>
308  (
309  "0",
310  phi.dimensions()/dimTime,
311  Zero
312  )
313  );
314 }
315 
316 
317 template<class Type>
319 (
320  const VolField<Type>& vf
321 )
322 {
324  (
325  "meshPhi",
326  mesh(),
328  );
329 }
330 
331 
332 template<class Type>
334 (
335  const VolField<Type>&,
336  const label patchi
337 )
338 {
339  return tmp<scalarField>
340  (
341  new scalarField(mesh().boundary()[patchi].size(), 0)
342  );
343 }
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace fv
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace Foam
353 
354 // ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:310
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:118
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:115
A class for managing temporary objects.
Definition: tmp.H:55
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
const dimensionSet dimVol
const dimensionSet dimVolume
const dimensionSet dimArea
faceListList boundary(nPatches)
labelList fv(nPoints)