basicThermo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "basicThermo.H"
30 #include "stringOps.H"
31 #include "wordIOList.H"
36 #include "fixedJumpFvPatchFields.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(basicThermo, 0);
46  defineRunTimeSelectionTable(basicThermo, fvMesh);
47  defineRunTimeSelectionTable(basicThermo, fvMeshDictPhase);
48 }
49 
50 const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
51 
52 const Foam::wordList Foam::basicThermo::componentHeader4
53 ({
54  "type",
55  "mixture",
56  "properties",
57  "energy"
58 });
59 
60 const Foam::wordList Foam::basicThermo::componentHeader7
61 ({
62  "type",
63  "mixture",
64  "transport",
65  "thermo",
66  "equationOfState",
67  "specie",
68  "energy"
69 });
70 
71 
72 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
73 
75 (
76  Ostream& os,
77  const wordList& cmptNames,
78  const wordList& thermoNames
79 )
80 {
81  const int nCmpt = cmptNames.size();
82 
83  // Build a table of constituent parts by split name into constituent parts
84  // - remove incompatible entries from the list
85  // - note: row-0 contains the names of constituent parts (ie, the header)
86 
87  DynamicList<wordList> outputTbl;
88  outputTbl.resize(thermoNames.size()+1);
89 
90  label rowi = 0;
91 
92  // Header
93  outputTbl[rowi] = cmptNames;
94  if (!outputTbl[rowi].empty())
95  {
96  ++rowi;
97  }
98 
99  for (const word& thermoName : thermoNames)
100  {
101  outputTbl[rowi] = basicThermo::splitThermoName(thermoName, nCmpt);
102  if (!outputTbl[rowi].empty())
103  {
104  ++rowi;
105  }
106  }
107 
108  if (rowi > 1)
109  {
110  outputTbl.resize(rowi);
111  Foam::printTable(outputTbl, os);
112  }
113 
114  return os;
115 }
116 
117 
118 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
119 
120 Foam::word Foam::basicThermo::makeThermoName
121 (
122  const dictionary& thermoTypeDict,
123  const wordList*& cmptHeaderPtr
124 )
125 {
126  if (thermoTypeDict.found("properties"))
127  {
128  if (cmptHeaderPtr)
129  {
130  cmptHeaderPtr = &(componentHeader4);
131  }
132 
133  return word
134  (
135  thermoTypeDict.get<word>("type") + '<'
136  + thermoTypeDict.get<word>("mixture") + '<'
137  + thermoTypeDict.get<word>("properties") + ','
138  + thermoTypeDict.get<word>("energy") + ">>"
139  );
140  }
141  else
142  {
143  if (cmptHeaderPtr)
144  {
145  cmptHeaderPtr = &(componentHeader7);
146  }
147 
148  return word
149  (
150  thermoTypeDict.get<word>("type") + '<'
151  + thermoTypeDict.get<word>("mixture") + '<'
152  + thermoTypeDict.get<word>("transport") + '<'
153  + thermoTypeDict.get<word>("thermo") + '<'
154  + thermoTypeDict.get<word>("equationOfState") + '<'
155  + thermoTypeDict.get<word>("specie") + ">>,"
156  + thermoTypeDict.get<word>("energy") + ">>>"
157  );
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
163 
165 {
166  const volScalarField::Boundary& tbf = this->T_.boundaryField();
167 
168  wordList hbt(tbf.size());
169 
170  forAll(tbf, patchi)
171  {
172  if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
173  {
174  const auto& pf =
175  dynamic_cast<const fixedJumpFvPatchScalarField&>
176  (
177  tbf[patchi]
178  );
179 
180  hbt[patchi] = pf.interfaceFieldType();
181  }
182  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
183  {
184  const auto& pf =
185  dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
186  (
187  tbf[patchi]
188  );
189 
190  hbt[patchi] = pf.interfaceFieldType();
191  }
192  }
193 
194  return hbt;
195 }
196 
197 
199 {
200  const volScalarField::Boundary& tbf = this->T_.boundaryField();
201 
202  wordList hbt(tbf.types());
203 
204  forAll(tbf, patchi)
205  {
206  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
207  {
208  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
209  }
210  else if
211  (
212  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
213  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
214  )
215  {
216  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
217  }
218  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
219  {
220  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
221  }
222  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
223  {
225  }
226  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
227  {
229  }
230  }
231 
232  return hbt;
233 }
234 
235 
236 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
237 
238 Foam::volScalarField& Foam::basicThermo::lookupOrConstruct
239 (
240  const fvMesh& mesh,
241  const word& fieldName,
242  bool& isOwner
243 )
244 {
245  auto* ptr = mesh.objectRegistry::getObjectPtr<volScalarField>(fieldName);
246 
247  isOwner = !ptr;
248 
249  if (!ptr)
250  {
251  ptr = new volScalarField
252  (
253  IOobject
254  (
255  fieldName,
256  mesh.time().timeName(),
257  mesh,
260  ),
261  mesh
262  );
263 
264  // Transfer ownership of this object to the objectRegistry
265  ptr->store();
266  }
267 
268  return *ptr;
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
273 
275 (
276  const fvMesh& mesh,
277  const word& phaseName
278 )
279 :
281  (
282  IOobject
283  (
284  phasePropertyName(dictName, phaseName),
285  mesh.time().constant(),
286  mesh,
287  IOobject::MUST_READ_IF_MODIFIED,
288  IOobject::NO_WRITE
289  )
290  ),
291 
292  phaseName_(phaseName),
293 
294  pOwner_(false),
295  TOwner_(false),
296  dpdt_(getOrDefault<bool>("dpdt", true)),
297 
298  p_(lookupOrConstruct(mesh, "p", pOwner_)),
299  T_(lookupOrConstruct(mesh, phasePropertyName("T"), TOwner_)),
300 
301  alpha_
302  (
303  IOobject
304  (
305  phasePropertyName("thermo:alpha"),
306  mesh.time().timeName(),
307  mesh,
308  IOobject::READ_IF_PRESENT,
309  IOobject::NO_WRITE
310  ),
311  mesh,
312  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
313  )
314 {
315  this->readIfPresent("updateT", TOwner_); // Manual override
316 }
317 
318 
320 (
321  const fvMesh& mesh,
322  const dictionary& dict,
323  const word& phaseName
324 )
325 :
327  (
328  IOobject
329  (
330  phasePropertyName(dictName, phaseName),
331  mesh.time().constant(),
332  mesh,
333  IOobject::NO_READ,
334  IOobject::NO_WRITE
335  ),
336  dict
337  ),
338 
339  phaseName_(phaseName),
340 
341  pOwner_(false),
342  TOwner_(false),
343  dpdt_(getOrDefault<bool>("dpdt", true)),
344 
345  p_(lookupOrConstruct(mesh, "p", pOwner_)),
346  T_(lookupOrConstruct(mesh, phasePropertyName("T"), TOwner_)),
347 
348  alpha_
349  (
350  IOobject
351  (
352  phasePropertyName("thermo:alpha"),
353  mesh.time().timeName(),
354  mesh,
355  IOobject::NO_READ,
356  IOobject::NO_WRITE
357  ),
358  mesh,
359  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
360  )
361 {
362  this->readIfPresent("updateT", TOwner_); // Manual override
363 }
364 
365 
367 (
368  const fvMesh& mesh,
369  const word& phaseName,
370  const word& dictionaryName
371 )
372 :
374  (
375  IOobject
376  (
377  dictionaryName,
378  mesh.time().constant(),
379  mesh,
380  IOobject::MUST_READ_IF_MODIFIED,
381  IOobject::NO_WRITE
382  )
383  ),
384 
385  phaseName_(phaseName),
386 
387  pOwner_(false),
388  TOwner_(false),
389  dpdt_(getOrDefault<bool>("dpdt", true)),
390 
391  p_(lookupOrConstruct(mesh, "p", pOwner_)),
392  T_(lookupOrConstruct(mesh, "T", TOwner_)),
393 
394  alpha_
395  (
396  IOobject
397  (
398  phasePropertyName("thermo:alpha"),
399  mesh.time().timeName(),
400  mesh,
401  IOobject::READ_IF_PRESENT,
402  IOobject::NO_WRITE
403  ),
404  mesh,
405  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
406  )
407 {
408  this->readIfPresent("updateT", TOwner_); // Manual override
409 
410  if (debug)
411  {
412  Pout<< "Constructed shared thermo : mesh:" << mesh.name()
413  << " phase:" << phaseName
414  << " dictionary:" << dictionaryName
415  << " T:" << T_.name()
416  << " updateT:" << TOwner_
417  << " alphaName:" << alpha_.name()
418  << endl;
419  }
420 }
421 
422 
423 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
424 
426 (
427  const fvMesh& mesh,
428  const word& phaseName
429 )
430 {
431  return New<basicThermo>(mesh, phaseName);
432 }
433 
434 
435 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
436 
438 {
439  if (pOwner_)
440  {
441  db().checkOut(p_.name());
442  }
443 
444  if (TOwner_)
445  {
446  db().checkOut(T_.name());
447  }
448 }
449 
450 
451 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
452 
454 (
455  const fvPatchScalarField& pf
456 )
457 {
458  const basicThermo* thermo = pf.db().findObject<basicThermo>(dictName);
459 
460  if (thermo)
461  {
462  return *thermo;
463  }
464 
465  HashTable<const basicThermo*> thermos =
466  pf.db().lookupClass<basicThermo>();
467 
469  {
470  thermo = iter.val();
471  if
472  (
473  &(thermo->he().internalField())
474  == &(pf.internalField())
475  )
476  {
477  return *thermo;
478  }
479  }
480 
481  return pf.db().lookupObject<basicThermo>(dictName);
482 }
483 
484 
486 (
487  const string& app,
488  const word& a
489 ) const
490 {
491  if (!(he().name() == phasePropertyName(a)))
492  {
494  << "Supported energy type is " << phasePropertyName(a)
495  << ", thermodynamics package provides " << he().name()
496  << exit(FatalError);
497  }
498 }
499 
501 (
502  const string& app,
503  const word& a,
504  const word& b
505 ) const
506 {
507  if
508  (
509  !(
510  he().name() == phasePropertyName(a)
511  || he().name() == phasePropertyName(b)
512  )
513  )
514  {
516  << "Supported energy types: " << phasePropertyName(a)
517  << " and " << phasePropertyName(b)
518  << ", thermodynamics package provides " << he().name()
519  << exit(FatalError);
520  }
521 }
522 
524 (
525  const string& app,
526  const word& a,
527  const word& b,
528  const word& c
529 ) const
530 {
531  if
532  (
533  !(
534  he().name() == phasePropertyName(a)
535  || he().name() == phasePropertyName(b)
536  || he().name() == phasePropertyName(c)
537  )
538  )
539  {
541  << "Supported energy types: " << phasePropertyName(a)
542  << ", " << phasePropertyName(b)
543  << " and " << phasePropertyName(c)
544  << ", thermodynamics package provides " << he().name()
545  << exit(FatalError);
546  }
547 }
548 
550 (
551  const string& app,
552  const word& a,
553  const word& b,
554  const word& c,
555  const word& d
556 ) const
557 {
558  if
559  (
560  !(
561  he().name() == phasePropertyName(a)
562  || he().name() == phasePropertyName(b)
563  || he().name() == phasePropertyName(c)
564  || he().name() == phasePropertyName(d)
565  )
566  )
567  {
569  << "Supported energy types: " << phasePropertyName(a)
570  << ", " << phasePropertyName(b)
571  << ", " << phasePropertyName(c)
572  << " and " << phasePropertyName(d)
573  << ", thermodynamics package provides " << he().name()
574  << exit(FatalError);
575  }
576 }
577 
578 
580 (
581  const std::string& thermoName,
582  const int nExpectedCmpts
583 )
584 {
585  // Split on ",<>" but include space for good measure.
586  // Splits things like
587  // "hePsiThermo<pureMixture<const<hConst<perfectGas<specie>>,enthalpy>>>"
588 
589  const auto parsed = stringOps::splitAny<std::string>(thermoName, " ,<>");
590  const int nParsed(parsed.size());
591 
592  wordList cmpts;
593 
594  if (!nExpectedCmpts || nParsed == nExpectedCmpts)
595  {
596  cmpts.resize(nParsed);
597 
598  auto iter = cmpts.begin();
599  for (const auto& sub : parsed)
600  {
601  *iter = word(sub.str());
602  ++iter;
603  }
604  }
605 
606  return cmpts;
607 }
608 
611 {
612  return p_;
613 }
614 
617 {
618  return p_;
619 }
620 
623 {
624  return T_;
625 }
626 
629 {
630  return T_;
631 }
632 
635 {
636  return alpha_;
637 }
638 
640 const Foam::scalarField& Foam::basicThermo::alpha(const label patchi) const
641 {
642  return alpha_.boundaryField()[patchi];
643 }
644 
645 
647 {
648  return regIOobject::read();
649 }
650 
651 
652 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:430
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:132
virtual bool read()
Read object.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
static Ostream & printThermoNames(Ostream &os, const wordList &cmptNames, const wordList &thermoNames)
Print (filtered) table of thermo names, splits on " ,<>".
Definition: basicThermo.C:68
Ostream & printTable(const UList< wordList > &tbl, List< std::string::size_type > &columnWidths, Ostream &os, bool headerSeparator=true)
Print a List of wordList as a table.
Definition: wordIOList.C:40
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
virtual word thermoName() const =0
Return the name of the thermo physics.
void resize(const label len)
Alter addressable list size, allocating new space if required while recovering old content...
Definition: DynamicListI.H:346
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:627
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:615
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
static autoPtr< Thermo > New(const fvMesh &, const word &phaseName=word::null)
Generic New for each of the related thermodynamics packages.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:603
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
psiReactionThermo & thermo
Definition: createFields.H:28
static const char *const typeName
Typename for Field.
Definition: Field.H:86
word timeName
Definition: getTimeIndex.H:3
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
wordList heBoundaryBaseTypes()
Return the enthalpy/internal energy field boundary base types by interrogating the temperature field ...
Definition: basicThermo.C:157
basicThermo(const basicThermo &)=delete
No copy construct.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
static wordList splitThermoName(const std::string &thermoName, const int nExpectedCmpts)
Split thermo package name into a list of components names.
Definition: basicThermo.C:573
fvPatchField< scalar > fvPatchScalarField
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
PtrList< solidThermo > thermos(solidRegions.size())
wordList heBoundaryTypes()
Return the enthalpy/internal energy field boundary types by interrogating the temperature field bound...
Definition: basicThermo.C:191
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:441
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
List< word > wordList
A List of words.
Definition: fileName.H:58
volScalarField alpha_
Laminar thermal diffusivity [kg/m/s].
Definition: basicThermo.H:166
bool TOwner_
Temperature created and stored by this instance.
Definition: basicThermo.H:143
const word & name() const
Return reference to name.
Definition: fvMesh.H:388
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const dimensionedScalar c
Speed of light in a vacuum.
Automatically write from objectRegistry::writeObject()
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:639
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Definition: basicThermo.C:447
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:479
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
volScalarField & T_
Temperature [K].
Definition: basicThermo.H:161
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157