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-2025 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 (
163 ) const
164 {
165  return Cv_()*T;
166 }
167 
168 
170 (
171  const scalarField& T,
172  const label patchi
173 ) const
174 {
175  return Cv_.boundaryField()[patchi]*T;
176 }
177 
178 
180 (
182  const fvSource& model,
183  const volScalarField::Internal& source
184 ) const
185 {
187  return tmp<volScalarField::Internal>(nullptr);
188 }
189 
190 
192 (
193  const scalarField& T,
194  const fvSource& model,
195  const scalarField& source,
196  const labelUList& cells
197 ) const
198 {
200  return tmp<scalarField>(nullptr);
201 }
202 
203 
205 {
207  return tmp<volScalarField>(nullptr);
208 }
209 
210 
212 (
213  const volScalarField& p,
214  const volScalarField& T
215 ) const
216 {
218  return tmp<volScalarField>(nullptr);
219 }
220 
221 
223 (
226 ) const
227 {
229  return tmp<volScalarField::Internal>(nullptr);
230 }
231 
232 
234 (
235  const scalarField& T,
236  const label patchi
237 ) const
238 {
240  return tmp<scalarField>(nullptr);
241 }
242 
243 
245 (
246  const scalarField& T,
247  const labelList& cells
248 ) const
249 {
251  return tmp<scalarField>(nullptr);
252 }
253 
254 
256 {
258  return tmp<volScalarField>(nullptr);
259 }
260 
261 
263 (
264  const volScalarField& p,
265  const volScalarField& T
266 ) const
267 {
269  return tmp<volScalarField>(nullptr);
270 }
271 
272 
274 (
277 ) const
278 {
280  return tmp<volScalarField::Internal>(nullptr);
281 }
282 
283 
285 (
286  const scalarField& T,
287  const label patchi
288 ) const
289 {
291  return tmp<scalarField>(nullptr);
292 }
293 
294 
296 (
297  const scalarField& T,
298  const labelList& cells
299 ) const
300 {
302  return tmp<scalarField>(nullptr);
303 }
304 
305 
307 (
308  const scalarField& T,
309  const label patchi
310 ) const
311 {
312  return Cv_.boundaryField()[patchi];
313 }
314 
315 
317 (
318  const scalarField& T,
319  const label patchi
320 ) const
321 {
322  return Cv_.boundaryField()[patchi];
323 }
324 
325 
327 (
328  const scalarField& T,
329  const label patchi
330 ) const
331 {
332  return Cv_.boundaryField()[patchi];
333 }
334 
335 
337 (
338  const volScalarField& h,
339  const volScalarField& p,
340  const volScalarField& T0
341 ) const
342 {
344  return tmp<volScalarField>(nullptr);
345 }
346 
347 
349 (
350  const scalarField& he,
351  const scalarField& T0,
352  const labelList& cells
353 ) const
354 {
356  return tmp<scalarField>(nullptr);
357 }
358 
359 
361 (
362  const scalarField& he,
363  const scalarField& T0,
364  const label patchi
365 ) const
366 {
368  return tmp<scalarField>(nullptr);
369 }
370 
371 
373 {
375  return volVectorField::null();
376 }
377 
378 
380 {
381  T_ = e_/Cv_;
382 }
383 
384 
385 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< Type, GeoMesh, PrimitiveField > & 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:504
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:96
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
thermo he()
volScalarField & p
scalar T0
Definition: createFields.H:22