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-2019 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 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(basicThermo, 0);
44  defineRunTimeSelectionTable(basicThermo, fvMesh);
45 }
46 
47 const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 {
54  const volScalarField::Boundary& tbf =
55  this->T_.boundaryField();
56 
57  wordList hbt(tbf.size(), word::null);
58 
59  forAll(tbf, patchi)
60  {
61  if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
62  {
63  const fixedJumpFvPatchScalarField& pf =
64  dynamic_cast<const fixedJumpFvPatchScalarField&>(tbf[patchi]);
65 
66  hbt[patchi] = pf.interfaceFieldType();
67  }
68  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
69  {
70  const fixedJumpAMIFvPatchScalarField& pf =
71  dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
72  (
73  tbf[patchi]
74  );
75 
76  hbt[patchi] = pf.interfaceFieldType();
77  }
78  }
79 
80  return hbt;
81 }
82 
83 
85 {
86  const volScalarField::Boundary& tbf =
87  this->T_.boundaryField();
88 
89  wordList hbt = tbf.types();
90 
91  forAll(tbf, patchi)
92  {
93  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
94  {
95  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
96  }
97  else if
98  (
99  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
100  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
101  || isA<gradientEnergyCalculatedTemperatureFvPatchScalarField>
102  (
103  tbf[patchi]
104  )
105  )
106  {
107  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
108  }
109  else if
110  (
111  isA<mixedFvPatchScalarField>(tbf[patchi])
112  || isA<mixedEnergyCalculatedTemperatureFvPatchScalarField>
113  (
114  tbf[patchi]
115  )
116  )
117  {
118  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
119  }
120  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
121  {
123  }
124  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
125  {
127  }
128  }
129 
130  return hbt;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
137 (
138  const fvMesh& mesh,
139  const char* name
140 ) const
141 {
142  if (!mesh.objectRegistry::foundObject<volScalarField>(name))
143  {
144  volScalarField* fPtr
145  (
146  new volScalarField
147  (
148  IOobject
149  (
150  name,
151  mesh.time().timeName(),
152  mesh,
155  ),
156  mesh
157  )
158  );
159 
160  // Transfer ownership of this object to the objectRegistry
161  fPtr->store(fPtr);
162  }
163 
164  return mesh.objectRegistry::lookupObjectRef<volScalarField>(name);
165 }
166 
167 
169 (
170  const fvMesh& mesh,
171  const word& phaseName
172 )
173 :
175  (
176  IOobject
177  (
178  phasePropertyName(dictName, phaseName),
179  mesh.time().constant(),
180  mesh,
183  )
184  ),
185 
186  phaseName_(phaseName),
187 
188  T_
189  (
190  IOobject
191  (
192  phasePropertyName("T"),
193  mesh.time().timeName(),
194  mesh,
197  ),
198  mesh
199  ),
200 
201  alpha_
202  (
203  IOobject
204  (
205  phasePropertyName("thermo:alpha"),
206  mesh.time().timeName(),
207  mesh,
210  ),
211  mesh,
212  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
213  ),
214 
215  dpdt_(lookupOrDefault<Switch>("dpdt", true))
216 {}
217 
218 
220 (
221  const fvMesh& mesh,
222  const dictionary& dict,
223  const word& phaseName
224 )
225 :
227  (
228  IOobject
229  (
230  phasePropertyName(dictName, phaseName),
231  mesh.time().constant(),
232  mesh,
235  ),
236  dict
237  ),
238 
239  phaseName_(phaseName),
240 
241  T_
242  (
243  IOobject
244  (
245  phasePropertyName("T"),
246  mesh.time().timeName(),
247  mesh,
250  ),
251  mesh
252  ),
253 
254  alpha_
255  (
256  IOobject
257  (
258  phasePropertyName("thermo:alpha"),
259  mesh.time().timeName(),
260  mesh,
263  ),
264  mesh,
265  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
266  )
267 {}
268 
269 
270 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
271 
273 (
274  const fvMesh& mesh,
275  const word& phaseName
276 )
277 {
278  return New<basicThermo>(mesh, phaseName);
279 }
280 
281 
282 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
283 
285 {}
286 
287 
288 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289 
291 (
292  const fvPatchScalarField& pf
293 )
294 {
295  if (pf.db().foundObject<basicThermo>(dictName))
296  {
297  return pf.db().lookupObject<basicThermo>(dictName);
298  }
299  else
300  {
302  pf.db().lookupClass<basicThermo>();
303 
304  for
305  (
307  iter != thermos.end();
308  ++iter
309  )
310  {
311  if
312  (
313  &(iter()->he().internalField())
314  == &(pf.internalField())
315  )
316  {
317  return *iter();
318  }
319  }
320  }
321 
322  return pf.db().lookupObject<basicThermo>(dictName);
323 }
324 
325 
327 (
328  const string& app,
329  const word& a
330 ) const
331 {
332  if (!(he().name() == phasePropertyName(a)))
333  {
335  << "Supported energy type is " << phasePropertyName(a)
336  << ", thermodynamics package provides " << he().name()
337  << exit(FatalError);
338  }
339 }
340 
342 (
343  const string& app,
344  const word& a,
345  const word& b
346 ) const
347 {
348  if
349  (
350  !(
351  he().name() == phasePropertyName(a)
352  || he().name() == phasePropertyName(b)
353  )
354  )
355  {
357  << "Supported energy types are " << phasePropertyName(a)
358  << " and " << phasePropertyName(b)
359  << ", thermodynamics package provides " << he().name()
360  << exit(FatalError);
361  }
362 }
363 
365 (
366  const string& app,
367  const word& a,
368  const word& b,
369  const word& c
370 ) const
371 {
372  if
373  (
374  !(
375  he().name() == phasePropertyName(a)
376  || he().name() == phasePropertyName(b)
377  || he().name() == phasePropertyName(c)
378  )
379  )
380  {
382  << "Supported energy types are " << phasePropertyName(a)
383  << ", " << phasePropertyName(b)
384  << " and " << phasePropertyName(c)
385  << ", thermodynamics package provides " << he().name()
386  << exit(FatalError);
387  }
388 }
389 
391 (
392  const string& app,
393  const word& a,
394  const word& b,
395  const word& c,
396  const word& d
397 ) const
398 {
399  if
400  (
401  !(
402  he().name() == phasePropertyName(a)
403  || he().name() == phasePropertyName(b)
404  || he().name() == phasePropertyName(c)
405  || he().name() == phasePropertyName(d)
406  )
407  )
408  {
410  << "Supported energy types are " << phasePropertyName(a)
411  << ", " << phasePropertyName(b)
412  << ", " << phasePropertyName(c)
413  << " and " << phasePropertyName(d)
414  << ", thermodynamics package provides " << he().name()
415  << exit(FatalError);
416  }
417 }
418 
419 
421 (
422  const word& thermoName,
423  const int nCmpt
424 )
425 {
426  wordList cmpts(nCmpt);
427 
428  string::size_type beg=0, end=0, endb=0, endc=0;
429  int i = 0;
430 
431  while
432  (
433  (endb = thermoName.find('<', beg)) != string::npos
434  || (endc = thermoName.find(',', beg)) != string::npos
435  )
436  {
437  if (endb == string::npos)
438  {
439  end = endc;
440  }
441  else if ((endc = thermoName.find(',', beg)) != string::npos)
442  {
443  end = min(endb, endc);
444  }
445  else
446  {
447  end = endb;
448  }
449 
450  if (beg < end)
451  {
452  cmpts[i] = thermoName.substr(beg, end-beg);
453  cmpts[i++].replaceAll(">","");
454 
455  // If the number of number of components in the name
456  // is greater than nCmpt return an empty list
457  if (i == nCmpt)
458  {
459  return wordList();
460  }
461  }
462  beg = end + 1;
463  }
464 
465  // If the number of number of components in the name is not equal to nCmpt
466  // return an empty list
467  if (i + 1 != nCmpt)
468  {
469  return wordList();
470  }
471 
472  if (beg < thermoName.size())
473  {
474  cmpts[i] = thermoName.substr(beg, string::npos);
475  cmpts[i].replaceAll(">","");
476  }
477 
478  return cmpts;
479 }
480 
481 
483 {
484  return T_;
485 }
486 
487 
489 {
490  return T_;
491 }
492 
493 
495 {
496  return alpha_;
497 }
498 
499 
501 {
502  return alpha_.boundaryField()[patchi];
503 }
504 
505 
507 {
508  return regIOobject::read();
509 }
510 
511 
512 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:284
#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
const word & name() const
Return name.
Definition: IOobject.H:303
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:494
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
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:482
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
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.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
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:52
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:421
A class for handling words, derived from string.
Definition: word.H:59
PtrList< solidThermo > thermos(solidRegions.size())
wordList heBoundaryTypes()
Return the enthalpy/internal energy field boundary types.
Definition: basicThermo.C:84
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static const word null
An empty word.
Definition: word.H:77
const word dictName("particleTrackDict")
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
volScalarField & lookupOrConstruct(const fvMesh &mesh, const char *name) const
Definition: basicThermo.C:137
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
static Table::iterator lookupThermo(const dictionary &thermoTypeDict, Table *tablePtr, const int nCmpt, const char *cmptNames[], const word &thermoTypeName)
Generic lookup for thermodynamics package thermoTypeName.
defineTypeNameAndDebug(combustionModel, 0)
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.
basicThermo(const basicThermo &)
Construct as copy (not implemented)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:167
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:506
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:327
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.
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332