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-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 "basicThermo.H"
33 #include "fixedJumpFvPatchFields.H"
35 
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 fvPatchScalarField& pf
83 )
84 {
85  if (pf.db().foundObject<basicThermo>(physicalProperties::typeName))
86  {
87  return pf.db().lookupObject<basicThermo>(physicalProperties::typeName);
88  }
89  else
90  {
92  pf.db().lookupClass<basicThermo>();
93 
94  for
95  (
97  iter != thermos.end();
98  ++iter
99  )
100  {
101  if
102  (
103  &(iter()->he().internalField())
104  == &(pf.internalField())
105  )
106  {
107  return *iter();
108  }
109  }
110  }
111 
112  return pf.db().lookupObject<basicThermo>(physicalProperties::typeName);
113 }
114 
115 
117 (
118  const word& thermoName,
119  const int nCmpt
120 )
121 {
122  wordList cmpts(nCmpt);
123 
124  string::size_type beg=0, end=0, endb=0, endc=0;
125  int i = 0;
126 
127  while
128  (
129  (endb = thermoName.find('<', beg)) != string::npos
130  || (endc = thermoName.find(',', beg)) != string::npos
131  )
132  {
133  if (endb == string::npos)
134  {
135  end = endc;
136  }
137  else if ((endc = thermoName.find(',', beg)) != string::npos)
138  {
139  end = min(endb, endc);
140  }
141  else
142  {
143  end = endb;
144  }
145 
146  if (beg < end)
147  {
148  cmpts[i] = thermoName.substr(beg, end-beg);
149  cmpts[i++].replaceAll(">","");
150 
151  // If the number of number of components in the name
152  // is greater than nCmpt return an empty list
153  if (i == nCmpt)
154  {
155  return wordList();
156  }
157  }
158  beg = end + 1;
159  }
160 
161  // If the number of number of components in the name is not equal to nCmpt
162  // return an empty list
163  if (i + 1 != nCmpt)
164  {
165  return wordList();
166  }
167 
168  if (beg < thermoName.size())
169  {
170  cmpts[i] = thermoName.substr(beg, string::npos);
171  cmpts[i].replaceAll(">","");
172  }
173 
174  return cmpts;
175 }
176 
177 
179 (
180  const word& thermoName
181 )
182 {
183  const wordList components(splitThermoName(thermoName, 5));
184 
185  return List<Pair<word>>
186  {
187  {"transport", components[0]},
188  {"thermo", components[1]},
189  {"equationOfState", components[2]},
190  {"specie", components[3]},
191  {"energy", components[4]}
192  };
193 }
194 
195 
196 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
197 
199 {
200  const volScalarField::Boundary& tbf = T().boundaryField();
201 
202  wordList hbt(tbf.size(), word::null);
203 
204  forAll(tbf, patchi)
205  {
206  if (tbf[patchi].overridesConstraint())
207  {
208  hbt[patchi] = tbf[patchi].patch().type();
209  }
210  }
211 
212  return hbt;
213 }
214 
215 
217 {
218  const volScalarField::Boundary& tbf = T().boundaryField();
219 
220  wordList hbt = tbf.types();
221 
222  forAll(tbf, patchi)
223  {
224  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
225  {
226  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
227  }
228  else if
229  (
230  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
231  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
232  || isA<gradientEnergyCalculatedTemperatureFvPatchScalarField>
233  (
234  tbf[patchi]
235  )
236  )
237  {
238  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
239  }
240  else if
241  (
242  isA<mixedFvPatchScalarField>(tbf[patchi])
243  || isA<mixedEnergyCalculatedTemperatureFvPatchScalarField>
244  (
245  tbf[patchi]
246  )
247  )
248  {
249  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
250  }
251  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
252  {
254  }
255  }
256 
257  return hbt;
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
262 
264 (
265  const fvMesh& mesh,
266  const word& phaseName
267 )
268 :
269  physicalProperties(mesh, phaseName),
270 
271  mesh_(mesh),
272 
273  phaseName_(phaseName),
274 
275  T_
276  (
277  IOobject
278  (
279  phasePropertyName("T", phaseName),
280  mesh.time().name(),
281  mesh,
282  IOobject::MUST_READ,
283  IOobject::AUTO_WRITE
284  ),
285  mesh
286  ),
287 
288  kappa_
289  (
290  IOobject
291  (
292  phasePropertyName("kappa", phaseName),
293  mesh.time().name(),
294  mesh,
295  IOobject::READ_IF_PRESENT,
296  IOobject::NO_WRITE
297  ),
298  mesh,
300  ),
301 
302  dpdt_(lookupOrDefault<Switch>("dpdt", true))
303 {}
304 
305 
306 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
307 
309 (
310  const fvMesh& mesh,
311  const word& phaseName
312 )
313 {
314  return New<basicThermo>(mesh, phaseName);
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
319 
321 {}
322 
323 
325 {}
326 
327 
328 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
329 
331 (
332  const string& app,
333  const word& a
334 ) const
335 {
336  if (!(he().name() == phasePropertyName(a)))
337  {
339  << "Supported energy type is " << phasePropertyName(a)
340  << ", thermodynamics package provides " << he().name()
341  << exit(FatalError);
342  }
343 }
344 
346 (
347  const string& app,
348  const word& a,
349  const word& b
350 ) const
351 {
352  if
353  (
354  !(
355  he().name() == phasePropertyName(a)
356  || he().name() == phasePropertyName(b)
357  )
358  )
359  {
361  << "Supported energy types are " << phasePropertyName(a)
362  << " and " << phasePropertyName(b)
363  << ", thermodynamics package provides " << he().name()
364  << exit(FatalError);
365  }
366 }
367 
368 
370 {
371  return T_;
372 }
373 
374 
376 {
377  return T_;
378 }
379 
380 
382 {
383  return kappa_;
384 }
385 
386 
388 {
389  return regIOobject::read();
390 }
391 
392 
393 // ************************************************************************* //
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
static const char *const typeName
Definition: Field.H:105
Generic GeometricBoundaryField class.
wordList types() const
Return a list of the patch field types.
Generic GeometricField class.
An STL-conforming iterator.
Definition: HashTable.H:429
An STL-conforming hash table.
Definition: HashTable.H:127
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
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:369
implementation(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: basicThermo.C:264
virtual const volScalarField & kappa() const
Thermal conductivity of mixture [W/m/K].
Definition: basicThermo.C:381
virtual ~implementation()
Destructor.
Definition: basicThermo.C:324
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:387
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Lookup the thermo associated with the given patch field.
Definition: basicThermo.C:81
static List< Pair< word > > thermoNameComponents(const word &thermoName)
Split name of thermo package into a list of named components names.
Definition: basicThermo.C:179
wordList heBoundaryTypes()
Enthalpy/internal energy field boundary types.
Definition: basicThermo.C:216
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
static word phasePropertyName(const word &name, const word &phaseName)
Name of a property for a given phase.
Definition: basicThermo.H:150
wordList heBoundaryBaseTypes()
Enthalpy/internal energy field boundary base types.
Definition: basicThermo.C:198
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:117
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const word & phaseName() const =0
Phase name.
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:331
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:320
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.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:361
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:145
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
bool foundObject(const word &name) const
Is the named Type in registry.
An abstract base class for physical properties.
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
virtual bool read()
Read object.
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:306
label patchi
volScalarField & b
Definition: createFields.H:27
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)
const dimensionSet dimEnergy
const dimensionSet dimLength
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
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)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112