liquidFilmThermo.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) 2013-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "liquidFilmThermo.H"
29 #include "demandDrivenData.H"
30 #include "thermoSingleLayer.H"
31 #include "SLGThermo.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace regionModels
40 {
41 namespace surfaceFilmModels
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(liquidFilmThermo, 0);
47 
49 (
50  filmThermoModel,
51  liquidFilmThermo,
52  dictionary
53 );
54 
55 
56 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
57 
59 {
60  if (!isA<thermoSingleLayer>(filmModel_))
61  {
63  << "Thermo model requires a " << thermoSingleLayer::typeName
64  << " film to supply the pressure and temperature, but "
65  << filmModel_.type() << " film model selected. "
66  << "Use the 'useReferenceValues' flag to employ reference "
67  << "pressure and temperature" << exit(FatalError);
68  }
69 
70  return refCast<const thermoSingleLayer>(filmModel_);
71 }
72 
73 
75 {
76  if (liquidPtr_ != nullptr)
77  {
78  return;
79  }
80 
81  dict.readEntry("liquid", name_);
82 
83  const SLGThermo* thermoPtr =
85 
86  if (thermoPtr)
87  {
88  // Retrieve from film thermo
89  ownLiquid_ = false;
90 
91  const SLGThermo& thermo = *thermoPtr;
92 
93  const label id = thermo.liquidId(name_);
94 
95  liquidPtr_ = &thermo.liquids().properties()[id];
96  }
97  else
98  {
99  // New liquid create
100  ownLiquid_ = true;
101 
102  liquidPtr_ =
103  liquidProperties::New(dict.optionalSubDict(name_ + "Coeffs")).ptr();
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
113  const dictionary& dict
114 )
115 :
116  filmThermoModel(typeName, film, dict),
117  name_("unknown_liquid"),
118  liquidPtr_(nullptr),
119  ownLiquid_(false),
120  useReferenceValues_(coeffDict_.get<bool>("useReferenceValues")),
121  pRef_(0.0),
122  TRef_(0.0)
123 {
124  initLiquid(coeffDict_);
125 
126  if (useReferenceValues_)
127  {
128  coeffDict_.readEntry("pRef", pRef_);
129  coeffDict_.readEntry("TRef", TRef_);
130  }
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
137 {
138  if (ownLiquid_)
139  {
141  }
142 }
143 
144 
145 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
146 
148 {
149  return name_;
150 }
151 
152 
154 (
155  const scalar p,
156  const scalar T
157 ) const
158 {
159  return liquidPtr_->rho(p, T);
160 }
161 
162 
164 (
165  const scalar p,
166  const scalar T
167 ) const
168 {
169  return liquidPtr_->mu(p, T);
170 }
171 
172 
174 (
175  const scalar p,
176  const scalar T
177 ) const
178 {
179  return liquidPtr_->sigma(p, T);
180 }
181 
182 
184 (
185  const scalar p,
186  const scalar T
187 ) const
188 {
189  return liquidPtr_->Cp(p, T);
190 }
191 
192 
194 (
195  const scalar p,
196  const scalar T
197 ) const
198 {
199  return liquidPtr_->kappa(p, T);
200 }
201 
202 
203 scalar liquidFilmThermo::D
204 (
205  const scalar p,
206  const scalar T
207 ) const
208 {
209  return liquidPtr_->D(p, T);
210 }
211 
212 
214 (
215  const scalar p,
216  const scalar T
217 ) const
218 {
219  return liquidPtr_->hl(p, T);
220 }
221 
222 
224 (
225  const scalar p,
226  const scalar T
227 ) const
228 {
229  return liquidPtr_->pv(p, T);
230 }
231 
232 
233 scalar liquidFilmThermo::W() const
234 {
235  return liquidPtr_->W();
236 }
237 
238 
239 scalar liquidFilmThermo::Tb(const scalar p) const
240 {
241  return liquidPtr_->pvInvert(p);
242 }
243 
244 
246 {
248  (
249  new volScalarField
250  (
251  IOobject
252  (
253  type() + ":rho",
254  film().time().timeName(),
255  film().regionMesh(),
258  ),
259  film().regionMesh(),
261  extrapolatedCalculatedFvPatchScalarField::typeName
262  )
263  );
264 
265  scalarField& rho = trho.ref().primitiveFieldRef();
266 
268  {
269  rho = this->rho(pRef_, TRef_);
270  }
271  else
272  {
273  const thermoSingleLayer& film = thermoFilm();
274 
275  const volScalarField& T = film.T();
276  const volScalarField& p = film.pPrimary();
277 
278  forAll(rho, celli)
279  {
280  rho[celli] = this->rho(p[celli], T[celli]);
281  }
282  }
283 
284  trho.ref().correctBoundaryConditions();
285 
286  return trho;
287 }
288 
289 
291 {
293  (
294  new volScalarField
295  (
296  IOobject
297  (
298  type() + ":mu",
299  film().time().timeName(),
300  film().regionMesh(),
303  ),
304  film().regionMesh(),
306  extrapolatedCalculatedFvPatchScalarField::typeName
307  )
308  );
309 
310  scalarField& mu = tmu.ref().primitiveFieldRef();
311 
313  {
314  mu = this->mu(pRef_, TRef_);
315  }
316  else
317  {
318  const thermoSingleLayer& film = thermoFilm();
319 
320  const volScalarField& T = film.T();
321  const volScalarField& p = film.pPrimary();
322 
323  forAll(mu, celli)
324  {
325  mu[celli] = this->mu(p[celli], T[celli]);
326  }
327  }
328 
329  tmu.ref().correctBoundaryConditions();
330 
331  return tmu;
332 }
333 
334 
336 {
337  tmp<volScalarField> tsigma
338  (
339  new volScalarField
340  (
341  IOobject
342  (
343  type() + ":sigma",
344  film().time().timeName(),
345  film().regionMesh(),
348  ),
349  film().regionMesh(),
351  extrapolatedCalculatedFvPatchScalarField::typeName
352  )
353  );
354 
355  scalarField& sigma = tsigma.ref().primitiveFieldRef();
356 
358  {
359  sigma = this->sigma(pRef_, TRef_);
360  }
361  else
362  {
363  const thermoSingleLayer& film = thermoFilm();
364 
365  const volScalarField& T = film.T();
366  const volScalarField& p = film.pPrimary();
367 
368  forAll(sigma, celli)
369  {
370  sigma[celli] = this->sigma(p[celli], T[celli]);
371  }
372  }
373 
374  tsigma.ref().correctBoundaryConditions();
375 
376  return tsigma;
377 }
378 
379 
381 {
383  (
384  new volScalarField
385  (
386  IOobject
387  (
388  type() + ":Cp",
389  film().time().timeName(),
390  film().regionMesh(),
393  ),
394  film().regionMesh(),
396  extrapolatedCalculatedFvPatchScalarField::typeName
397  )
398  );
399 
400  scalarField& Cp = tCp.ref().primitiveFieldRef();
401 
403  {
404  Cp = this->Cp(pRef_, TRef_);
405  }
406  else
407  {
408  const thermoSingleLayer& film = thermoFilm();
409 
410  const volScalarField& T = film.T();
411  const volScalarField& p = film.pPrimary();
412 
413  forAll(Cp, celli)
414  {
415  Cp[celli] = this->Cp(p[celli], T[celli]);
416  }
417  }
418 
419  tCp.ref().correctBoundaryConditions();
420 
421  return tCp;
422 }
423 
424 
426 {
427  tmp<volScalarField> tkappa
428  (
429  new volScalarField
430  (
431  IOobject
432  (
433  type() + ":kappa",
434  film().time().timeName(),
435  film().regionMesh(),
438  ),
439  film().regionMesh(),
441  extrapolatedCalculatedFvPatchScalarField::typeName
442  )
443  );
444 
445  scalarField& kappa = tkappa.ref().primitiveFieldRef();
446 
448  {
449  kappa = this->kappa(pRef_, TRef_);
450  }
451  else
452  {
453  const thermoSingleLayer& film = thermoFilm();
454 
455  const volScalarField& T = film.T();
456  const volScalarField& p = film.pPrimary();
457 
458  forAll(kappa, celli)
459  {
460  kappa[celli] = this->kappa(p[celli], T[celli]);
461  }
462  }
463 
464  tkappa.ref().correctBoundaryConditions();
465 
466  return tkappa;
467 }
468 
469 
470 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
471 
472 } // End namespace surfaceFilmModels
473 } // End namespace regionModels
474 } // End namespace Foam
475 
476 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::thermoFilm
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
Definition: liquidFilmThermo.C:58
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::pRef_
scalar pRef_
Reference pressure [pa].
Definition: liquidFilmThermo.H:79
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::liquidProperties::New
static autoPtr< liquidProperties > New(const word &name)
Return a pointer to a new liquidProperties created from name.
Definition: liquidProperties.C:91
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::name
virtual const word & name() const
Return the specie name.
Definition: liquidFilmThermo.C:147
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::kappa
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: liquidFilmThermo.C:425
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::D
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
Definition: liquidFilmThermo.C:204
tCp
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::surfaceFilmModels::filmSubModelBase::film
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Definition: filmSubModelBaseI.H:39
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::sigma
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
Definition: liquidFilmThermo.C:335
SLGThermo.H
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::rho
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
Definition: liquidFilmThermo.C:245
Foam::thermophysicalProperties::W
scalar W() const
Molecular weight [kg/kmol].
Definition: thermophysicalPropertiesI.H:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
thermoSingleLayer.H
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::ownLiquid_
bool ownLiquid_
Flag to indicate that model owns the liquid object.
Definition: liquidFilmThermo.H:73
Foam::regionModels::regionModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:34
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
liquidFilmThermo.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Thermodynamic form of single-cell layer surface film model.
Definition: thermoSingleLayer.H:67
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::initLiquid
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
Definition: liquidFilmThermo.C:74
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::Cp
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: liquidFilmThermo.C:380
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::subModelBase::dict
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:113
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::T
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::mu
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
Definition: liquidFilmThermo.C:290
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::liquidPtr_
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
Definition: liquidFilmThermo.H:70
Foam::liquidProperties::pvInvert
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
Definition: liquidProperties.C:181
Foam::dimPower
const dimensionSet dimPower
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::hl
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
Definition: liquidFilmThermo.C:214
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::Tb
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
Definition: liquidFilmThermo.C:239
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::useReferenceValues_
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
Definition: liquidFilmThermo.H:76
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::regionModels::surfaceFilmModels::filmThermoModel
Base class for film thermo models.
Definition: filmThermoModel.H:56
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::TRef_
scalar TRef_
Reference temperature [K].
Definition: liquidFilmThermo.H:82
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::dictionary::optionalSubDict
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::W
virtual scalar W() const
Return molecular weight [kg/kmol].
Definition: liquidFilmThermo.C:233
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::name_
word name_
Liquid name.
Definition: liquidFilmThermo.H:67
Foam::regionModels::surfaceFilmModels::filmSubModelBase::filmModel_
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Definition: filmSubModelBase.H:65
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::~liquidFilmThermo
virtual ~liquidFilmThermo()
Destructor.
Definition: liquidFilmThermo.C:136
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::pv
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
Definition: liquidFilmThermo.C:224
extrapolatedCalculatedFvPatchFields.H
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::liquidFilmThermo
liquidFilmThermo(const liquidFilmThermo &)=delete
No copy construct.