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 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 :
48  solidThermo::composite(mesh, phaseName),
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  )
64 {
65  rho_ = readProperty<scalar>("rho", rho_.dimensions());
66 
67  if (readKappa)
68  {
69  kappa_ = readProperty<scalar>("kappa", kappa_.dimensions());
70  }
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  const fvMesh& mesh,
79  const word& phaseName
80 )
81 :
82  constSolidThermo(mesh, true, phaseName)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  return Cv_;
97 }
98 
99 
101 {
102  return e_;
103 }
104 
105 
107 {
108  return e_;
109 }
110 
111 
113 (
114  const scalarField& T,
115  const labelList& cells
116 ) const
117 {
118  return scalarField(Cv_, cells)*T;
119 }
120 
121 
123 (
124  const volScalarField& p,
125  const volScalarField& T
126 ) const
127 {
128  return Cv_*T;
129 }
130 
131 
133 (
134  const scalarField& T,
135  const label patchi
136 ) const
137 {
138  return Cv_.boundaryField()[patchi]*T;
139 }
140 
141 
143 {
145  return tmp<volScalarField>(nullptr);
146 }
147 
148 
150 (
151  const volScalarField& p,
152  const volScalarField& T
153 ) const
154 {
156  return tmp<volScalarField>(nullptr);
157 }
158 
159 
161 (
162  const scalarField& T,
163  const label patchi
164 ) const
165 {
167  return tmp<scalarField>(nullptr);
168 }
169 
170 
172 (
173  const scalarField& T,
174  const labelList& cells
175 ) const
176 {
178  return tmp<scalarField>(nullptr);
179 }
180 
181 
183 {
185  return tmp<volScalarField>(nullptr);
186 }
187 
188 
190 (
191  const volScalarField& p,
192  const volScalarField& T
193 ) const
194 {
196  return tmp<volScalarField>(nullptr);
197 }
198 
199 
201 (
202  const scalarField& T,
203  const label patchi
204 ) const
205 {
207  return tmp<scalarField>(nullptr);
208 }
209 
210 
212 (
213  const scalarField& T,
214  const labelList& cells
215 ) const
216 {
218  return tmp<scalarField>(nullptr);
219 }
220 
221 
223 {
225  return tmp<volScalarField>(nullptr);
226 }
227 
228 
230 (
231  const volScalarField& h,
232  const volScalarField& p,
233  const volScalarField& T0
234 ) const
235 {
237  return tmp<volScalarField>(nullptr);
238 }
239 
240 
242 (
243  const scalarField& he,
244  const scalarField& T0,
245  const labelList& cells
246 ) const
247 {
249  return tmp<scalarField>(nullptr);
250 }
251 
252 
254 (
255  const scalarField& he,
256  const scalarField& T0,
257  const label patchi
258 ) const
259 {
261  return tmp<scalarField>(nullptr);
262 }
263 
264 
266 {
267  return Cv_;
268 }
269 
270 
272 {
273  return Cv_;
274 }
275 
276 
278 (
279  const scalarField& T,
280  const label patchi
281 ) const
282 {
283  return Cv_.boundaryField()[patchi];
284 }
285 
286 
288 (
289  const scalarField& T,
290  const label patchi
291 ) const
292 {
293  return Cv_.boundaryField()[patchi];
294 }
295 
296 
298 (
299  const scalarField& T,
300  const label patchi
301 ) const
302 {
303  return Cv_.boundaryField()[patchi];
304 }
305 
306 
308 {
310  return volVectorField::null();
311 }
312 
313 
315 {
316  T_ = e_/Cv_;
317 }
318 
319 
321 {
322  return regIOobject::read();
323 }
324 
325 
326 // ************************************************************************* //
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
volScalarField kappa_
Thermal conductivity [W/m/K].
Definition: basicThermo.H:440
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
Uniform or non-uniform constant solid thermodynamic properties.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual void correct()
Update properties.
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
virtual ~constSolidThermo()
Destructor.
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
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 tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
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 bool read()
Read thermophysicalProperties dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
virtual bool read()
Read object.
volScalarField rho_
Density field [kg/m^3].
Definition: solidThermo.H:139
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:56
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:353
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
const dimensionSet dimEnergy
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
thermo he()
volScalarField & p
scalar T0
Definition: createFields.H:22