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-2021 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"
37 
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(basicThermo, 0);
46  defineRunTimeSelectionTable(basicThermo, fvMesh);
47 }
48 
49 
50 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51 
53 (
54  const fvMesh& mesh,
55  const char* name
56 )
57 {
58  if (!mesh.objectRegistry::foundObject<volScalarField>(name))
59  {
60  volScalarField* fPtr
61  (
62  new volScalarField
63  (
64  IOobject
65  (
66  name,
67  mesh.time().timeName(),
68  mesh,
71  ),
72  mesh
73  )
74  );
75 
76  // Transfer ownership of this object to the objectRegistry
77  fPtr->store(fPtr);
78  }
79 
80  return mesh.objectRegistry::lookupObjectRef<volScalarField>(name);
81 }
82 
83 
85 (
86  const fvPatchScalarField& pf
87 )
88 {
90  {
91  return pf.db().lookupObject<basicThermo>(dictName);
92  }
93  else
94  {
96  pf.db().lookupClass<basicThermo>();
97 
98  for
99  (
101  iter != thermos.end();
102  ++iter
103  )
104  {
105  if
106  (
107  &(iter()->he().internalField())
108  == &(pf.internalField())
109  )
110  {
111  return *iter();
112  }
113  }
114  }
115 
116  return pf.db().lookupObject<basicThermo>(dictName);
117 }
118 
119 
121 (
122  const word& thermoName,
123  const int nCmpt
124 )
125 {
126  wordList cmpts(nCmpt);
127 
128  string::size_type beg=0, end=0, endb=0, endc=0;
129  int i = 0;
130 
131  while
132  (
133  (endb = thermoName.find('<', beg)) != string::npos
134  || (endc = thermoName.find(',', beg)) != string::npos
135  )
136  {
137  if (endb == string::npos)
138  {
139  end = endc;
140  }
141  else if ((endc = thermoName.find(',', beg)) != string::npos)
142  {
143  end = min(endb, endc);
144  }
145  else
146  {
147  end = endb;
148  }
149 
150  if (beg < end)
151  {
152  cmpts[i] = thermoName.substr(beg, end-beg);
153  cmpts[i++].replaceAll(">","");
154 
155  // If the number of number of components in the name
156  // is greater than nCmpt return an empty list
157  if (i == nCmpt)
158  {
159  return wordList();
160  }
161  }
162  beg = end + 1;
163  }
164 
165  // If the number of number of components in the name is not equal to nCmpt
166  // return an empty list
167  if (i + 1 != nCmpt)
168  {
169  return wordList();
170  }
171 
172  if (beg < thermoName.size())
173  {
174  cmpts[i] = thermoName.substr(beg, string::npos);
175  cmpts[i].replaceAll(">","");
176  }
177 
178  return cmpts;
179 }
180 
181 
183 (
184  const word& thermoName
185 )
186 {
187  const wordList components(splitThermoName(thermoName, 5));
188 
189  return List<Pair<word>>
190  {
191  {"transport", components[0]},
192  {"thermo", components[1]},
193  {"equationOfState", components[2]},
194  {"specie", components[3]},
195  {"energy", components[4]}
196  };
197 }
198 
199 
200 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
201 
203 {
204  const volScalarField::Boundary& tbf = T().boundaryField();
205 
206  wordList hbt(tbf.size(), word::null);
207 
208  forAll(tbf, patchi)
209  {
210  if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
211  {
212  const fixedJumpFvPatchScalarField& pf =
213  dynamic_cast<const fixedJumpFvPatchScalarField&>(tbf[patchi]);
214 
215  hbt[patchi] = pf.interfaceFieldType();
216  }
217  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
218  {
219  const fixedJumpAMIFvPatchScalarField& pf =
220  dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
221  (
222  tbf[patchi]
223  );
224 
225  hbt[patchi] = pf.interfaceFieldType();
226  }
227  }
228 
229  return hbt;
230 }
231 
232 
234 {
235  const volScalarField::Boundary& tbf = T().boundaryField();
236 
237  wordList hbt = tbf.types();
238 
239  forAll(tbf, patchi)
240  {
241  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
242  {
243  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
244  }
245  else if
246  (
247  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
248  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
249  || isA<gradientEnergyCalculatedTemperatureFvPatchScalarField>
250  (
251  tbf[patchi]
252  )
253  )
254  {
255  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
256  }
257  else if
258  (
259  isA<mixedFvPatchScalarField>(tbf[patchi])
260  || isA<mixedEnergyCalculatedTemperatureFvPatchScalarField>
261  (
262  tbf[patchi]
263  )
264  )
265  {
266  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
267  }
268  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
269  {
271  }
272  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
273  {
275  }
276  }
277 
278  return hbt;
279 }
280 
281 
282 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
283 
285 (
286  const fvMesh& mesh,
287  const word& phaseName
288 )
289 :
291  (
292  IOobject
293  (
294  phasePropertyName(dictName, phaseName),
295  mesh.time().constant(),
296  mesh,
299  )
300  ),
301 
302  phaseName_(phaseName),
303 
304  T_
305  (
306  IOobject
307  (
308  phasePropertyName("T", phaseName),
309  mesh.time().timeName(),
310  mesh,
313  ),
314  mesh
315  ),
316 
317  alpha_
318  (
319  IOobject
320  (
321  phasePropertyName("thermo:alpha", phaseName),
322  mesh.time().timeName(),
323  mesh,
326  ),
327  mesh,
328  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
329  ),
330 
331  dpdt_(lookupOrDefault<Switch>("dpdt", true))
332 {}
333 
334 
336 (
337  const fvMesh& mesh,
338  const dictionary& dict,
339  const word& phaseName
340 )
341 :
343  (
344  IOobject
345  (
346  phasePropertyName(dictName, phaseName),
347  mesh.time().constant(),
348  mesh,
351  ),
352  dict
353  ),
354 
355  phaseName_(phaseName),
356 
357  T_
358  (
359  IOobject
360  (
361  phasePropertyName("T", phaseName),
362  mesh.time().timeName(),
363  mesh,
366  ),
367  mesh
368  ),
369 
370  alpha_
371  (
372  IOobject
373  (
374  phasePropertyName("thermo:alpha", phaseName),
375  mesh.time().timeName(),
376  mesh,
379  ),
380  mesh,
381  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
382  )
383 {}
384 
385 
386 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
387 
389 (
390  const fvMesh& mesh,
391  const word& phaseName
392 )
393 {
394  return New<basicThermo>(mesh, phaseName);
395 }
396 
397 
398 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
399 
401 {}
402 
403 
405 {}
406 
407 
408 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
409 
411 (
412  const string& app,
413  const word& a
414 ) const
415 {
416  if (!(he().name() == phasePropertyName(a)))
417  {
419  << "Supported energy type is " << phasePropertyName(a)
420  << ", thermodynamics package provides " << he().name()
421  << exit(FatalError);
422  }
423 }
424 
426 (
427  const string& app,
428  const word& a,
429  const word& b
430 ) const
431 {
432  if
433  (
434  !(
435  he().name() == phasePropertyName(a)
436  || he().name() == phasePropertyName(b)
437  )
438  )
439  {
441  << "Supported energy types are " << phasePropertyName(a)
442  << " and " << phasePropertyName(b)
443  << ", thermodynamics package provides " << he().name()
444  << exit(FatalError);
445  }
446 }
447 
448 
450 {
451  return T_;
452 }
453 
454 
456 {
457  return T_;
458 }
459 
460 
462 {
463  return alpha_;
464 }
465 
466 
468 (
469  const label patchi
470 ) const
471 {
472  return alpha_.boundaryField()[patchi];
473 }
474 
475 
477 {
478  return regIOobject::read();
479 }
480 
481 
482 // ************************************************************************* //
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:400
implementation(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: basicThermo.C:285
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
virtual bool read()
Read object.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static const char *const typeName
Definition: Field.H:105
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual const volScalarField & alpha() const
Return the thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:461
bool foundObject(const word &name) const
Is the named Type found?
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
static autoPtr< Thermo > New(const fvMesh &, const word &phaseName=word::null)
Generic New for each of the related thermodynamics packages.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual ~implementation()
Destructor.
Definition: basicThermo.C:404
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
An STL-conforming iterator.
Definition: HashTable.H:426
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
wordList heBoundaryBaseTypes()
Return the enthalpy/internal energy field boundary base types.
Definition: basicThermo.C:202
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:121
A class for handling words, derived from string.
Definition: word.H:59
wordList heBoundaryTypes()
Return the enthalpy/internal energy field boundary types.
Definition: basicThermo.C:233
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:449
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
static const word null
An empty word.
Definition: word.H:77
static const zero Zero
Definition: zero.H:97
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
An STL-conforming hash table.
Definition: HashTable.H:61
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
thermo he()
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:147
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Lookup the thermo associated with the given patch field.
Definition: basicThermo.C:85
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:476
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:411
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
static List< Pair< word > > thermoNameComponents(const word &thermoName)
Split name of thermo package into a list of named components names.
Definition: basicThermo.C:183
const word dictName("noiseDict")
Namespace for OpenFOAM.
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:53
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:343