constSolidThermo.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) 2022-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 "constSolidThermo.H"
28 
29 /* * * * * * * * * * * * * * * Private Static Data * * * * * * * * * * * * * */
30 
31 namespace Foam
32 {
36 }
37 
38 
39 // * * * * * * * * * * * * * * Protected Constructors * * * * * * * * * * * //
40 
42 (
43  const fvMesh& mesh,
44  const bool readKappa,
45  const word& phaseName
46 )
47 :
49  Cv_(readProperty<scalar>("Cv", dimEnergy/dimMass/dimTemperature)),
50  e_
51  (
52  IOobject
53  (
54  phasePropertyName("e"),
55  mesh.time().name(),
56  mesh,
57  IOobject::NO_READ,
58  IOobject::NO_WRITE
59  ),
60  Cv_*T_,
61  this->heBoundaryTypes(),
62  this->heBoundaryBaseTypes(),
63  this->heSourcesTypes()
64  )
65 {
66  rho_ = readProperty<scalar>("rho", rho_.dimensions());
67 
68  if (readKappa)
69  {
70  kappa_ = readProperty<scalar>("kappa", kappa_.dimensions());
71  }
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 (
79  const fvMesh& mesh,
80  const word& phaseName
81 )
82 :
83  constSolidThermo(mesh, true, phaseName)
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
98  return tmp<volScalarField>(nullptr);
99 }
100 
101 
103 {
105  return tmp<scalarField>(nullptr);
106 }
107 
108 
110 {
111  return e_;
112 }
113 
114 
116 {
117  return e_;
118 }
119 
120 
122 {
123  return Cv_;
124 }
125 
126 
128 {
129  return Cv_;
130 }
131 
132 
134 {
135  return Cv_;
136 }
137 
138 
140 (
141  const scalarField& T,
142  const labelList& cells
143 ) const
144 {
145  return scalarField(Cv_, cells)*T;
146 }
147 
148 
150 (
151  const volScalarField& p,
152  const volScalarField& T
153 ) const
154 {
155  return Cv_*T;
156 }
157 
158 
160 (
161  const scalarField& T,
162  const label patchi
163 ) const
164 {
165  return Cv_.boundaryField()[patchi]*T;
166 }
167 
168 
170 (
171  const scalarField& T,
172  const fvSource& source
173 ) const
174 {
176  return tmp<scalarField>(nullptr);
177 }
178 
179 
181 {
183  return tmp<volScalarField>(nullptr);
184 }
185 
186 
188 (
189  const volScalarField& p,
190  const volScalarField& T
191 ) const
192 {
194  return tmp<volScalarField>(nullptr);
195 }
196 
197 
199 (
200  const scalarField& T,
201  const label patchi
202 ) const
203 {
205  return tmp<scalarField>(nullptr);
206 }
207 
208 
210 (
211  const scalarField& T,
212  const labelList& cells
213 ) const
214 {
216  return tmp<scalarField>(nullptr);
217 }
218 
219 
221 {
223  return tmp<volScalarField>(nullptr);
224 }
225 
226 
228 (
229  const volScalarField& p,
230  const volScalarField& T
231 ) const
232 {
234  return tmp<volScalarField>(nullptr);
235 }
236 
237 
239 (
240  const scalarField& T,
241  const label patchi
242 ) const
243 {
245  return tmp<scalarField>(nullptr);
246 }
247 
248 
250 (
251  const scalarField& T,
252  const labelList& cells
253 ) const
254 {
256  return tmp<scalarField>(nullptr);
257 }
258 
259 
261 (
262  const scalarField& T,
263  const label patchi
264 ) const
265 {
266  return Cv_.boundaryField()[patchi];
267 }
268 
269 
271 (
272  const scalarField& T,
273  const label patchi
274 ) const
275 {
276  return Cv_.boundaryField()[patchi];
277 }
278 
279 
281 (
282  const scalarField& T,
283  const label patchi
284 ) const
285 {
286  return Cv_.boundaryField()[patchi];
287 }
288 
289 
291 (
292  const volScalarField& h,
293  const volScalarField& p,
294  const volScalarField& T0
295 ) const
296 {
298  return tmp<volScalarField>(nullptr);
299 }
300 
301 
303 (
304  const scalarField& he,
305  const scalarField& T0,
306  const labelList& cells
307 ) const
308 {
310  return tmp<scalarField>(nullptr);
311 }
312 
313 
315 (
316  const scalarField& he,
317  const scalarField& T0,
318  const label patchi
319 ) const
320 {
322  return tmp<scalarField>(nullptr);
323 }
324 
325 
327 {
329  return volVectorField::null();
330 }
331 
332 
334 {
335  T_ = e_/Cv_;
336 }
337 
338 
339 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Wrapper around a thermo which also constructs the physical properties dictionary.
volScalarField kappa_
Thermal conductivity [W/m/K].
Definition: basicThermo.H:470
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
Uniform or non-uniform constant solid thermodynamic properties.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual void correct()
Update properties.
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
virtual ~constSolidThermo()
Destructor.
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg].
constSolidThermo(const fvMesh &, const bool readKappa, const word &phaseName=word::null)
Construct from mesh and phase name.
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
virtual const volVectorField & Kappa() const
Anisotropic thermal conductivity [W/m/K].
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Base class for finite volume sources.
Definition: fvSource.H:52
volScalarField rho_
Density field [kg/m^3].
Definition: rhoThermo.H:95
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
const cellShapeList & cells
const dimensionedScalar h
Planck constant.
Namespace for OpenFOAM.
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
const dimensionSet dimEnergy
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
thermo he()
volScalarField & p
scalar T0
Definition: createFields.H:22