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-2018 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"
31 #include "fixedJumpFvPatchFields.H"
35 
36 
37 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(basicThermo, 0);
42  defineRunTimeSelectionTable(basicThermo, fvMesh);
43 }
44 
45 const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  const volScalarField::Boundary& tbf =
53  this->T_.boundaryField();
54 
55  wordList hbt(tbf.size(), word::null);
56 
57  forAll(tbf, patchi)
58  {
59  if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
60  {
61  const fixedJumpFvPatchScalarField& pf =
62  dynamic_cast<const fixedJumpFvPatchScalarField&>(tbf[patchi]);
63 
64  hbt[patchi] = pf.interfaceFieldType();
65  }
66  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
67  {
68  const fixedJumpAMIFvPatchScalarField& pf =
69  dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
70  (
71  tbf[patchi]
72  );
73 
74  hbt[patchi] = pf.interfaceFieldType();
75  }
76  }
77 
78  return hbt;
79 }
80 
81 
83 {
84  const volScalarField::Boundary& tbf =
85  this->T_.boundaryField();
86 
87  wordList hbt = tbf.types();
88 
89  forAll(tbf, patchi)
90  {
91  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
92  {
93  hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
94  }
95  else if
96  (
97  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
98  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
99  )
100  {
101  hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
102  }
103  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
104  {
105  hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
106  }
107  else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
108  {
110  }
111  else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
112  {
114  }
115  else if (tbf[patchi].type() == "energyRegionCoupledFvPatchScalarField")
116  {
117  hbt[patchi] = "energyRegionCoupledFvPatchScalarField";
118  }
119  }
120 
121  return hbt;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 (
129  const fvMesh& mesh,
130  const char* name
131 ) const
132 {
133  if (!mesh.objectRegistry::foundObject<volScalarField>(name))
134  {
135  volScalarField* fPtr
136  (
137  new volScalarField
138  (
139  IOobject
140  (
141  name,
142  mesh.time().timeName(),
143  mesh,
146  ),
147  mesh
148  )
149  );
150 
151  // Transfer ownership of this object to the objectRegistry
152  fPtr->store(fPtr);
153  }
154 
155  return mesh.objectRegistry::lookupObjectRef<volScalarField>(name);
156 }
157 
158 
160 (
161  const fvMesh& mesh,
162  const word& phaseName
163 )
164 :
166  (
167  IOobject
168  (
169  phasePropertyName(dictName, phaseName),
170  mesh.time().constant(),
171  mesh,
174  )
175  ),
176 
177  phaseName_(phaseName),
178 
179  p_(lookupOrConstruct(mesh, "p")),
180 
181  T_
182  (
183  IOobject
184  (
185  phasePropertyName("T"),
186  mesh.time().timeName(),
187  mesh,
190  ),
191  mesh
192  ),
193 
194  alpha_
195  (
196  IOobject
197  (
198  phasePropertyName("thermo:alpha"),
199  mesh.time().timeName(),
200  mesh,
203  ),
204  mesh,
206  (
207  "zero",
208  dimensionSet(1, -1, -1, 0, 0),
209  Zero
210  )
211  ),
212 
213  dpdt_(lookupOrDefault<Switch>("dpdt", true))
214 {}
215 
216 
218 (
219  const fvMesh& mesh,
220  const dictionary& dict,
221  const word& phaseName
222 )
223 :
225  (
226  IOobject
227  (
228  phasePropertyName(dictName, phaseName),
229  mesh.time().constant(),
230  mesh,
233  ),
234  dict
235  ),
236 
237  phaseName_(phaseName),
238 
239  p_(lookupOrConstruct(mesh, "p")),
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,
266  (
267  "zero",
268  dimensionSet(1, -1, -1, 0, 0),
269  Zero
270  )
271  )
272 {}
273 
274 
275 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
276 
278 (
279  const fvMesh& mesh,
280  const word& phaseName
281 )
282 {
283  return New<basicThermo>(mesh, phaseName);
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
288 
290 {}
291 
292 
293 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
294 
296 (
297  const fvPatchScalarField& pf
298 )
299 {
300  if (pf.db().foundObject<basicThermo>(dictName))
301  {
302  return pf.db().lookupObject<basicThermo>(dictName);
303  }
304  else
305  {
307  pf.db().lookupClass<basicThermo>();
308 
309  for
310  (
312  iter != thermos.end();
313  ++iter
314  )
315  {
316  if
317  (
318  &(iter()->he().internalField())
319  == &(pf.internalField())
320  )
321  {
322  return *iter();
323  }
324  }
325  }
326 
327  return pf.db().lookupObject<basicThermo>(dictName);
328 }
329 
330 
332 (
333  const string& app,
334  const word& a
335 ) const
336 {
337  if (!(he().name() == phasePropertyName(a)))
338  {
340  << "Supported energy type is " << phasePropertyName(a)
341  << ", thermodynamics package provides " << he().name()
342  << exit(FatalError);
343  }
344 }
345 
347 (
348  const string& app,
349  const word& a,
350  const word& b
351 ) const
352 {
353  if
354  (
355  !(
356  he().name() == phasePropertyName(a)
357  || he().name() == phasePropertyName(b)
358  )
359  )
360  {
362  << "Supported energy types are " << phasePropertyName(a)
363  << " and " << phasePropertyName(b)
364  << ", thermodynamics package provides " << he().name()
365  << exit(FatalError);
366  }
367 }
368 
370 (
371  const string& app,
372  const word& a,
373  const word& b,
374  const word& c
375 ) const
376 {
377  if
378  (
379  !(
380  he().name() == phasePropertyName(a)
381  || he().name() == phasePropertyName(b)
382  || he().name() == phasePropertyName(c)
383  )
384  )
385  {
387  << "Supported energy types are " << phasePropertyName(a)
388  << ", " << phasePropertyName(b)
389  << " and " << phasePropertyName(c)
390  << ", thermodynamics package provides " << he().name()
391  << exit(FatalError);
392  }
393 }
394 
396 (
397  const string& app,
398  const word& a,
399  const word& b,
400  const word& c,
401  const word& d
402 ) const
403 {
404  if
405  (
406  !(
407  he().name() == phasePropertyName(a)
408  || he().name() == phasePropertyName(b)
409  || he().name() == phasePropertyName(c)
410  || he().name() == phasePropertyName(d)
411  )
412  )
413  {
415  << "Supported energy types are " << phasePropertyName(a)
416  << ", " << phasePropertyName(b)
417  << ", " << phasePropertyName(c)
418  << " and " << phasePropertyName(d)
419  << ", thermodynamics package provides " << he().name()
420  << exit(FatalError);
421  }
422 }
423 
424 
426 (
427  const word& thermoName,
428  const int nCmpt
429 )
430 {
431  wordList cmpts(nCmpt);
432 
433  string::size_type beg=0, end=0, endb=0, endc=0;
434  int i = 0;
435 
436  while
437  (
438  (endb = thermoName.find('<', beg)) != string::npos
439  || (endc = thermoName.find(',', beg)) != string::npos
440  )
441  {
442  if (endb == string::npos)
443  {
444  end = endc;
445  }
446  else if ((endc = thermoName.find(',', beg)) != string::npos)
447  {
448  end = min(endb, endc);
449  }
450  else
451  {
452  end = endb;
453  }
454 
455  if (beg < end)
456  {
457  cmpts[i] = thermoName.substr(beg, end-beg);
458  cmpts[i++].replaceAll(">","");
459 
460  // If the number of number of components in the name
461  // is greater than nCmpt return an empty list
462  if (i == nCmpt)
463  {
464  return wordList();
465  }
466  }
467  beg = end + 1;
468  }
469 
470  // If the number of number of components in the name is not equal to nCmpt
471  // return an empty list
472  if (i + 1 != nCmpt)
473  {
474  return wordList();
475  }
476 
477  if (beg < thermoName.size())
478  {
479  cmpts[i] = thermoName.substr(beg, string::npos);
480  cmpts[i].replaceAll(">","");
481  }
482 
483  return cmpts;
484 }
485 
486 
488 {
489  return p_;
490 }
491 
492 
494 {
495  return p_;
496 }
497 
498 
500 {
501  return T_;
502 }
503 
504 
506 {
507  return T_;
508 }
509 
510 
512 {
513  return alpha_;
514 }
515 
516 
518 {
519  return alpha_.boundaryField()[patchi];
520 }
521 
522 
524 {
525  return regIOobject::read();
526 }
527 
528 
529 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:51
virtual ~basicThermo()
Destructor.
Definition: basicThermo.C:289
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:297
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:110
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:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static const char *const typeName
Definition: Field.H:94
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:511
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:499
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
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:115
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:487
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:50
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Definition: basicThermo.C:426
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:82
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:419
An STL-conforming hash table.
Definition: HashTable.H:62
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:128
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
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
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:197
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:523
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:332
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:347