Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
48}
49
50const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
51
52const Foam::wordList Foam::basicThermo::componentHeader4
53({
54 "type",
55 "mixture",
56 "properties",
57 "energy"
58});
59
60const 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
120Foam::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();
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 (
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
238Foam::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 (
254 (
255 fieldName,
256 mesh.time().timeName(),
257 mesh,
261 ),
262 mesh
263 );
264
265 // Transfer ownership of this object to the objectRegistry
266 ptr->store();
267 }
269 return *ptr;
270}
271
272
273// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
274
276(
277 const fvMesh& mesh,
278 const word& phaseName
279)
280:
282 (
284 (
285 phasePropertyName(dictName, phaseName),
286 mesh.time().constant(),
287 mesh,
288 IOobject::READ_MODIFIED,
289 IOobject::NO_WRITE,
290 IOobject::REGISTER
291 )
292 ),
293
294 phaseName_(phaseName),
295
296 pOwner_(false),
297 TOwner_(false),
298 dpdt_(getOrDefault<bool>("dpdt", true)),
299
300 p_(lookupOrConstruct(mesh, "p", pOwner_)),
301 T_(lookupOrConstruct(mesh, phasePropertyName("T"), TOwner_)),
302
303 alpha_
304 (
306 (
307 phaseScopedName("thermo", "alpha"),
308 mesh.time().timeName(),
309 mesh,
310 IOobject::READ_IF_PRESENT,
311 IOobject::NO_WRITE,
312 IOobject::REGISTER
313 ),
314 mesh,
316 )
317{
318 this->readIfPresent("updateT", TOwner_); // Manual override
319}
320
321
323(
324 const fvMesh& mesh,
325 const dictionary& dict,
326 const word& phaseName
327)
328:
330 (
332 (
333 phasePropertyName(dictName, phaseName),
334 mesh.time().constant(),
335 mesh,
336 IOobject::NO_READ,
337 IOobject::NO_WRITE,
338 IOobject::REGISTER
339 ),
340 dict
341 ),
342
343 phaseName_(phaseName),
344
345 pOwner_(false),
346 TOwner_(false),
347 dpdt_(getOrDefault<bool>("dpdt", true)),
348
349 p_(lookupOrConstruct(mesh, "p", pOwner_)),
350 T_(lookupOrConstruct(mesh, phasePropertyName("T"), TOwner_)),
351
352 alpha_
353 (
355 (
356 phaseScopedName("thermo", "alpha"),
357 mesh.time().timeName(),
358 mesh,
359 IOobject::NO_READ,
360 IOobject::NO_WRITE,
361 IOobject::REGISTER
362 ),
363 mesh,
365 )
366{
367 this->readIfPresent("updateT", TOwner_); // Manual override
368}
369
370
372(
373 const fvMesh& mesh,
374 const word& phaseName,
375 const word& dictionaryName
376)
377:
379 (
381 (
382 dictionaryName,
383 mesh.time().constant(),
384 mesh,
385 IOobject::READ_MODIFIED,
386 IOobject::NO_WRITE,
387 IOobject::REGISTER
388 )
389 ),
390
391 phaseName_(phaseName),
392
393 pOwner_(false),
394 TOwner_(false),
395 dpdt_(getOrDefault<bool>("dpdt", true)),
396
397 p_(lookupOrConstruct(mesh, "p", pOwner_)),
398 T_(lookupOrConstruct(mesh, "T", TOwner_)),
399
400 alpha_
401 (
403 (
404 phaseScopedName("thermo", "alpha"),
405 mesh.time().timeName(),
406 mesh,
407 IOobject::READ_IF_PRESENT,
408 IOobject::NO_WRITE,
409 IOobject::REGISTER
410 ),
411 mesh,
412 dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
413 )
414{
415 this->readIfPresent("updateT", TOwner_); // Manual override
416
417 if (debug)
418 {
419 Pout<< "Constructed shared thermo : mesh:" << mesh.name()
420 << " phase:" << phaseName
421 << " dictionary:" << dictionaryName
422 << " T:" << T_.name()
423 << " updateT:" << TOwner_
424 << " alphaName:" << alpha_.name()
426 }
427}
428
429
430// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
431
433(
434 const fvMesh& mesh,
435 const word& phaseName
436)
438 return New<basicThermo>(mesh, phaseName);
439}
440
441
442// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
443
445{
446 if (pOwner_)
447 {
448 db().checkOut(p_.name());
449 }
450
451 if (TOwner_)
452 {
453 db().checkOut(T_.name());
454 }
455}
456
457
458// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
459
461(
462 const fvPatchScalarField& pf
463)
464{
465 const basicThermo* thermoPtr = pf.db().cfindObject<basicThermo>(dictName);
466
467 if (thermoPtr)
468 {
469 return *thermoPtr;
470 }
471
472 for (const basicThermo& thermo : pf.db().cobjects<basicThermo>())
473 {
474 if
475 (
476 &(thermo.he().internalField())
477 == &(pf.internalField())
478 )
479 {
480 return thermo;
481 }
483
484 // Failure
485 return pf.db().lookupObject<basicThermo>(dictName);
486}
487
488
490(
491 const string& app,
492 const word& a
493) const
494{
495 if (!(he().name() == phasePropertyName(a)))
496 {
498 << "Supported energy type is " << phasePropertyName(a)
499 << ", thermodynamics package provides " << he().name()
500 << exit(FatalError);
501 }
502}
503
505(
506 const string& app,
507 const word& a,
508 const word& b
509) const
510{
511 if
512 (
513 !(
514 he().name() == phasePropertyName(a)
515 || he().name() == phasePropertyName(b)
516 )
517 )
518 {
520 << "Supported energy types: " << phasePropertyName(a)
521 << " and " << phasePropertyName(b)
522 << ", thermodynamics package provides " << he().name()
523 << exit(FatalError);
524 }
525}
526
528(
529 const string& app,
530 const word& a,
531 const word& b,
532 const word& c
533) const
534{
535 if
536 (
537 !(
538 he().name() == phasePropertyName(a)
539 || he().name() == phasePropertyName(b)
540 || he().name() == phasePropertyName(c)
541 )
542 )
543 {
545 << "Supported energy types: " << phasePropertyName(a)
547 << " and " << phasePropertyName(c)
548 << ", thermodynamics package provides " << he().name()
549 << exit(FatalError);
550 }
551}
552
554(
555 const string& app,
556 const word& a,
557 const word& b,
558 const word& c,
559 const word& d
560) const
561{
562 if
563 (
564 !(
565 he().name() == phasePropertyName(a)
566 || he().name() == phasePropertyName(b)
567 || he().name() == phasePropertyName(c)
568 || he().name() == phasePropertyName(d)
569 )
570 )
571 {
573 << "Supported energy types: " << phasePropertyName(a)
574 << ", " << phasePropertyName(b)
575 << ", " << phasePropertyName(c)
576 << " and " << phasePropertyName(d)
577 << ", thermodynamics package provides " << he().name()
578 << exit(FatalError);
579 }
580}
581
582
584(
585 const std::string& thermoName,
586 const int nExpectedCmpts
587)
588{
589 // Split on ",<>" but include space for good measure.
590 // Splits things like
591 // "hePsiThermo<pureMixture<const<hConst<perfectGas<specie>>,enthalpy>>>"
592
593 const auto parsed = stringOps::splitAny(thermoName, " ,<>");
594 const int nParsed(parsed.size());
595
596 wordList cmpts;
597
598 if (!nExpectedCmpts || nParsed == nExpectedCmpts)
599 {
600 cmpts.resize(nParsed);
601
602 auto iter = cmpts.begin();
603 for (const auto& sub : parsed)
604 {
605 *iter = word(sub.str());
606 ++iter;
608 }
609
610 return cmpts;
611}
612
617}
618
621{
622 return p_;
623}
624
627{
628 return T_;
629}
630
635}
636
639{
640 return alpha_;
641}
642
644const Foam::scalarField& Foam::basicThermo::alpha(const label patchi) const
645{
646 return alpha_.boundaryField()[patchi];
647}
648
649
651{
652 return regIOobject::read();
653}
654
655
656// ************************************************************************* //
volScalarField & he
Definition YEEqn.H:52
bool empty() const noexcept
True if the list is empty.
Definition DLListBase.H:189
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void resize(const label len)
Alter addressable list size, allocating new space if required while recovering old content.
static const char *const typeName
Typename for Field.
Definition Field.H:93
wordList types() const
Return a list of the patch types.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
virtual word thermoName() const =0
Return the name of the thermo physics.
static Ostream & printThermoNames(Ostream &os, const wordList &cmptNames, const wordList &thermoNames)
Print (filtered) table of thermo names, splits on " ,<>".
Definition basicThermo.C:68
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
bool dpdt_
Include dpdt term in the enthalpy equation?
virtual const volScalarField & T() const
Temperature [K].
wordList heBoundaryTypes()
Return the enthalpy/internal energy field boundary types by interrogating the temperature field bound...
virtual volScalarField & p()
Pressure [Pa].
static word phasePropertyName(const word &name, const word &phaseName)
The phase property name as property.phase.
wordList heBoundaryBaseTypes()
Return the enthalpy/internal energy field boundary base types by interrogating the temperature field ...
static const word dictName
The dictionary name ("thermophysicalProperties").
bool pOwner_
Pressure created and stored by this instance.
bool TOwner_
Temperature created and stored by this instance.
volScalarField & p_
Pressure [Pa].
basicThermo(const basicThermo &)=delete
No copy construct.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
const word & phaseName_
Phase-name.
volScalarField & T_
Temperature [K].
word phaseScopedName(const std::string &scope, const word &name) const
The phase property name scoped as scope:property.phase.
static autoPtr< Thermo > New(const fvMesh &, const word &phaseName=word::null)
Generic New for each of the related thermodynamics packages.
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
virtual ~basicThermo()
Destructor.
static wordList splitThermoName(const std::string &thermoName, const int nExpectedCmpts)
Split thermo package name into a list of components names.
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
virtual bool read()
Read thermophysical properties dictionary.
volScalarField alpha_
Laminar thermal diffusivity [kg/m/s].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
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...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const objectRegistry & db() const
The associated objectRegistry.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
UPtrList< const Type > cobjects() const
Return unsorted list of objects with a class satisfying isA<Type> or isType<Type> (with Strict).
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...
virtual bool read()
Read object.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
auto & name
word timeName
Definition getTimeIndex.H:3
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
Foam::SubStrings splitAny(const std::string &str, const std::string &delim, std::string::size_type pos=0)
Split string into sub-strings using any characters in delimiter.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & b
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299