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"
37 
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(basicThermo, 0);
44  defineRunTimeSelectionTable(basicThermo, fvMesh);
45 }
46 
47 
48 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
49 
51 (
52  const fvMesh& mesh,
53  const char* name
54 )
55 {
56  if (!mesh.objectRegistry::foundObject<volScalarField>(name))
57  {
58  volScalarField* fPtr
59  (
60  new volScalarField
61  (
62  IOobject
63  (
64  name,
65  mesh.time().timeName(),
66  mesh,
69  ),
70  mesh
71  )
72  );
73 
74  // Transfer ownership of this object to the objectRegistry
75  fPtr->store(fPtr);
76  }
77 
78  return mesh.objectRegistry::lookupObjectRef<volScalarField>(name);
79 }
80 
81 
83 (
84  const fvPatchScalarField& pf
85 )
86 {
87  if (pf.db().foundObject<basicThermo>(physicalProperties::typeName))
88  {
89  return pf.db().lookupObject<basicThermo>(physicalProperties::typeName);
90  }
91  else
92  {
94  pf.db().lookupClass<basicThermo>();
95 
96  for
97  (
99  iter != thermos.end();
100  ++iter
101  )
102  {
103  if
104  (
105  &(iter()->he().internalField())
106  == &(pf.internalField())
107  )
108  {
109  return *iter();
110  }
111  }
112  }
113 
114  return pf.db().lookupObject<basicThermo>(physicalProperties::typeName);
115 }
116 
117 
119 (
120  const word& thermoName,
121  const int nCmpt
122 )
123 {
124  wordList cmpts(nCmpt);
125 
126  string::size_type beg=0, end=0, endb=0, endc=0;
127  int i = 0;
128 
129  while
130  (
131  (endb = thermoName.find('<', beg)) != string::npos
132  || (endc = thermoName.find(',', beg)) != string::npos
133  )
134  {
135  if (endb == string::npos)
136  {
137  end = endc;
138  }
139  else if ((endc = thermoName.find(',', beg)) != string::npos)
140  {
141  end = min(endb, endc);
142  }
143  else
144  {
145  end = endb;
146  }
147 
148  if (beg < end)
149  {
150  cmpts[i] = thermoName.substr(beg, end-beg);
151  cmpts[i++].replaceAll(">","");
152 
153  // If the number of number of components in the name
154  // is greater than nCmpt return an empty list
155  if (i == nCmpt)
156  {
157  return wordList();
158  }
159  }
160  beg = end + 1;
161  }
162 
163  // If the number of number of components in the name is not equal to nCmpt
164  // return an empty list
165  if (i + 1 != nCmpt)
166  {
167  return wordList();
168  }
169 
170  if (beg < thermoName.size())
171  {
172  cmpts[i] = thermoName.substr(beg, string::npos);
173  cmpts[i].replaceAll(">","");
174  }
175 
176  return cmpts;
177 }
178 
179 
181 (
182  const word& thermoName
183 )
184 {
185  const wordList components(splitThermoName(thermoName, 5));
186 
187  return List<Pair<word>>
188  {
189  {"transport", components[0]},
190  {"thermo", components[1]},
191  {"equationOfState", components[2]},
192  {"specie", components[3]},
193  {"energy", components[4]}
194  };
195 }
196 
197 
198 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
199 
201 {
202  const volScalarField::Boundary& tbf = T().boundaryField();
203 
204  wordList hbt(tbf.size(), word::null);
205 
206  forAll(tbf, patchi)
207  {
208  if (tbf[patchi].overridesConstraint())
209  {
210  hbt[patchi] = tbf[patchi].patch().type();
211  }
212  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
213  {
214  const fixedJumpAMIFvPatchScalarField& pf =
215  dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
216  (
217  tbf[patchi]
218  );
219 
220  hbt[patchi] = pf.interfaceFieldType();
221  }
222  }
223 
224  return hbt;
225 }
226 
227 
229 {
230  const volScalarField::Boundary& tbf = T().boundaryField();
231 
232  wordList hbt = tbf.types();
233 
234  forAll(tbf, patchi)
235  {
236  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
237  {
238  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
239  }
240  else if
241  (
242  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
243  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
244  || isA<gradientEnergyCalculatedTemperatureFvPatchScalarField>
245  (
246  tbf[patchi]
247  )
248  )
249  {
250  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
251  }
252  else if
253  (
254  isA<mixedFvPatchScalarField>(tbf[patchi])
255  || isA<mixedEnergyCalculatedTemperatureFvPatchScalarField>
256  (
257  tbf[patchi]
258  )
259  )
260  {
261  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
262  }
263  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
264  {
266  }
267  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
268  {
270  }
271  }
272 
273  return hbt;
274 }
275 
276 
277 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
278 
280 (
281  const fvMesh& mesh,
282  const word& phaseName
283 )
284 :
285  physicalProperties(mesh, phaseName),
286 
287  phaseName_(phaseName),
288 
289  T_
290  (
291  IOobject
292  (
293  phasePropertyName("T", phaseName),
294  mesh.time().timeName(),
295  mesh,
298  ),
299  mesh
300  ),
301 
302  kappa_
303  (
304  IOobject
305  (
306  phasePropertyName("thermo:kappa", phaseName),
307  mesh.time().timeName(),
308  mesh,
311  ),
312  mesh,
314  ),
315 
316  dpdt_(lookupOrDefault<Switch>("dpdt", true))
317 {}
318 
319 
320 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
321 
323 (
324  const fvMesh& mesh,
325  const word& phaseName
326 )
327 {
328  return New<basicThermo>(mesh, phaseName);
329 }
330 
331 
332 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
333 
335 {}
336 
337 
339 {}
340 
341 
342 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
343 
345 (
346  const string& app,
347  const word& a
348 ) const
349 {
350  if (!(he().name() == phasePropertyName(a)))
351  {
353  << "Supported energy type is " << phasePropertyName(a)
354  << ", thermodynamics package provides " << he().name()
355  << exit(FatalError);
356  }
357 }
358 
360 (
361  const string& app,
362  const word& a,
363  const word& b
364 ) const
365 {
366  if
367  (
368  !(
369  he().name() == phasePropertyName(a)
370  || he().name() == phasePropertyName(b)
371  )
372  )
373  {
375  << "Supported energy types are " << phasePropertyName(a)
376  << " and " << phasePropertyName(b)
377  << ", thermodynamics package provides " << he().name()
378  << exit(FatalError);
379  }
380 }
381 
382 
384 {
385  return T_;
386 }
387 
388 
390 {
391  return T_;
392 }
393 
394 
396 {
397  return kappa_;
398 }
399 
400 
402 {
403  return regIOobject::read();
404 }
405 
406 
407 // ************************************************************************* //
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:334
implementation(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: basicThermo.C:280
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
static const char *const typeName
Definition: Field.H:105
An abstract base class for physical properties.
wordList types() const
Return a list of the patch field types.
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
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
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:338
const dimensionSet dimLength
An STL-conforming iterator.
Definition: HashTable.H:426
const dimensionSet dimTime
wordList heBoundaryBaseTypes()
Enthalpy/internal energy field boundary base types.
Definition: basicThermo.C:200
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:119
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
wordList heBoundaryTypes()
Enthalpy/internal energy field boundary types.
Definition: basicThermo.C:228
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:383
static const word null
An empty word.
Definition: word.H:77
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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:40
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
thermo he()
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimEnergy
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
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:95
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:83
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:401
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:345
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
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:181
const dimensionSet dimTemperature
virtual const volScalarField & kappa() const
Thermal conductivity of mixture [W/m/K].
Definition: basicThermo.C:395
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:51
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:357