basicThermo.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-2024 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 "basicThermo.H"
33 #include "fixedJumpFvPatchFields.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
43 }
44 
45 
46 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
47 
49 (
50  const fvMesh& mesh,
51  const char* name
52 )
53 {
54  if (!mesh.objectRegistry::foundObject<volScalarField>(name))
55  {
56  volScalarField* fPtr
57  (
58  new volScalarField
59  (
60  IOobject
61  (
62  name,
63  mesh.time().name(),
64  mesh,
67  ),
68  mesh
69  )
70  );
71 
72  // Transfer ownership of this object to the objectRegistry
73  fPtr->store(fPtr);
74  }
75 
76  return mesh.objectRegistry::lookupObjectRef<volScalarField>(name);
77 }
78 
79 
81 (
82  const word& thermoName,
83  const int nCmpt
84 )
85 {
86  wordList cmpts(nCmpt);
87 
88  string::size_type beg=0, end=0, endb=0, endc=0;
89  int i = 0;
90 
91  while
92  (
93  (endb = thermoName.find('<', beg)) != string::npos
94  || (endc = thermoName.find(',', beg)) != string::npos
95  )
96  {
97  if (endb == string::npos)
98  {
99  end = endc;
100  }
101  else if ((endc = thermoName.find(',', beg)) != string::npos)
102  {
103  end = min(endb, endc);
104  }
105  else
106  {
107  end = endb;
108  }
109 
110  if (beg < end)
111  {
112  cmpts[i] = thermoName.substr(beg, end-beg);
113  cmpts[i++].replaceAll(">","");
114 
115  // If the number of number of components in the name
116  // is greater than nCmpt return an empty list
117  if (i == nCmpt)
118  {
119  return wordList();
120  }
121  }
122  beg = end + 1;
123  }
124 
125  // If the number of number of components in the name is not equal to nCmpt
126  // return an empty list
127  if (i + 1 != nCmpt)
128  {
129  return wordList();
130  }
131 
132  if (beg < thermoName.size())
133  {
134  cmpts[i] = thermoName.substr(beg, string::npos);
135  cmpts[i].replaceAll(">","");
136  }
137 
138  return cmpts;
139 }
140 
141 
143 (
144  const word& thermoName
145 )
146 {
147  const wordList components(splitThermoName(thermoName, 5));
148 
149  return List<Pair<word>>
150  {
151  {"transport", components[0]},
152  {"thermo", components[1]},
153  {"equationOfState", components[2]},
154  {"specie", components[3]},
155  {"energy", components[4]}
156  };
157 }
158 
159 
160 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
161 
163 {
164  const volScalarField::Boundary& tbf = T().boundaryField();
165 
166  wordList hbt(tbf.size(), word::null);
167 
168  forAll(tbf, patchi)
169  {
170  if (tbf[patchi].overridesConstraint())
171  {
172  hbt[patchi] = tbf[patchi].patch().type();
173  }
174  }
175 
176  return hbt;
177 }
178 
179 
181 {
182  const volScalarField::Boundary& tbf = T().boundaryField();
183 
184  wordList hbt = tbf.types();
185 
186  forAll(tbf, patchi)
187  {
188  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
189  {
190  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
191  }
192  else if
193  (
194  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
195  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
196  || isA<gradientEnergyCalculatedTemperatureFvPatchScalarField>
197  (
198  tbf[patchi]
199  )
200  )
201  {
202  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
203  }
204  else if
205  (
206  isA<mixedFvPatchScalarField>(tbf[patchi])
207  || isA<mixedEnergyCalculatedTemperatureFvPatchScalarField>
208  (
209  tbf[patchi]
210  )
211  )
212  {
213  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
214  }
215  else if (isA<jumpCyclicFvPatchScalarField>(tbf[patchi]))
216  {
217  hbt[patchi] = energyJumpFvPatchScalarField::typeName;
218  }
219  }
220 
221  return hbt;
222 }
223 
224 
226 {
227  const HashTable<word> tst = T().sources().types();
228 
229  HashTable<word> hst;
230  forAllConstIter(typename HashTable<word>, tst, iter)
231  {
232  hst.set(iter.key(), energyFvScalarFieldSource::typeName);
233  }
234 
235  return hst;
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
240 
242 (
243  const dictionary& dict,
244  const fvMesh& mesh,
245  const word& phaseName
246 )
247 :
248  mesh_(mesh),
249 
250  phaseName_(phaseName),
251 
252  T_
253  (
254  IOobject
255  (
256  phasePropertyName("T", phaseName),
257  mesh.time().name(),
258  mesh,
259  IOobject::MUST_READ,
260  IOobject::AUTO_WRITE
261  ),
262  mesh
263  ),
264 
265  kappa_
266  (
267  IOobject
268  (
269  phasePropertyName("kappa", phaseName),
270  mesh.time().name(),
271  mesh,
272  IOobject::READ_IF_PRESENT,
273  IOobject::NO_WRITE
274  ),
275  mesh,
277  ),
278 
279  dpdt_(dict.lookupOrDefault<Switch>("dpdt", true))
280 {}
281 
282 
283 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
284 
286 (
287  const fvMesh& mesh,
288  const word& phaseName
289 )
290 {
291  return New<basicThermo>(mesh, phaseName);
292 }
293 
294 
295 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
296 
298 {}
299 
300 
302 {}
303 
304 
305 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
306 
308 (
309  const string& app,
310  const word& a
311 ) const
312 {
313  if (!(he().name() == phasePropertyName(a)))
314  {
316  << "Supported energy type is " << phasePropertyName(a)
317  << ", thermodynamics package provides " << he().name()
318  << exit(FatalError);
319  }
320 }
321 
322 
324 (
325  const string& app,
326  const word& a,
327  const word& b
328 ) const
329 {
330  if
331  (
332  !(
333  he().name() == phasePropertyName(a)
334  || he().name() == phasePropertyName(b)
335  )
336  )
337  {
339  << "Supported energy types are " << phasePropertyName(a)
340  << " and " << phasePropertyName(b)
341  << ", thermodynamics package provides " << he().name()
342  << exit(FatalError);
343  }
344 }
345 
346 
348 {
349  return volScalarField::New(phasePropertyName("gamma"), Cp()/Cv());
350 }
351 
352 
354 (
355  const scalarField& T,
356  const label patchi
357 ) const
358 {
359  return Cp(T, patchi)/Cv(T, patchi);
360 }
361 
362 
364 {
365  return T_;
366 }
367 
368 
370 {
371  return T_;
372 }
373 
374 
376 {
377  return kappa_;
378 }
379 
380 
382 {}
383 
384 
385 // ************************************************************************* //
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricBoundaryField class.
wordList types() const
Return a list of the patch field types.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:363
virtual const volScalarField & kappa() const
Thermal conductivity of mixture [W/m/K].
Definition: basicThermo.C:375
virtual ~implementation()
Destructor.
Definition: basicThermo.C:301
implementation(const dictionary &, const fvMesh &, const word &)
Construct from dictionary, mesh and phase name.
Definition: basicThermo.C:242
virtual void read(const dictionary &)
Read thermophysical properties dictionary.
Definition: basicThermo.C:381
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
virtual const volScalarField & Cv() const =0
Heat capacity at constant volume [J/kg/K].
static List< Pair< word > > thermoNameComponents(const word &thermoName)
Split name of thermo package into a list of named components names.
Definition: basicThermo.C:143
wordList heBoundaryTypes()
Enthalpy/internal energy field boundary types.
Definition: basicThermo.C:180
HashTable< word > heSourcesTypes()
Enthalpy/internal energy field sources types.
Definition: basicThermo.C:225
static word phasePropertyName(const word &name, const word &phaseName)
Name of a property for a given phase.
Definition: basicThermo.H:154
wordList heBoundaryBaseTypes()
Enthalpy/internal energy field boundary base types.
Definition: basicThermo.C:162
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:81
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual const word & phaseName() const =0
Phase name.
tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Definition: basicThermo.C:347
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:308
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:297
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
static volScalarField & lookupOrConstruct(const fvMesh &mesh, const char *name)
Lookup and the named field, or construct it as MUST-READ if it is.
Definition: basicThermo.C:49
static autoPtr< Thermo > New(const fvMesh &, const word &phaseName=word::null)
Generic New for each of the related thermodynamics packages.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
const dimensionSet dimThermalConductivity
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict