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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "liquidFilmThermo.H"
30 #include "demandDrivenData.H"
31 #include "thermoSingleLayer.H"
32 #include "SLGThermo.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace regionModels
41 {
42 namespace surfaceFilmModels
43 {
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
47 defineTypeNameAndDebug(liquidFilmThermo, 0);
48 
50 (
51  filmThermoModel,
54 );
55 
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
60 {
61  if (!isA<thermoSingleLayer>(filmModel_))
62  {
64  << "Thermo model requires a " << thermoSingleLayer::typeName
65  << " film to supply the pressure and temperature, but "
66  << filmModel_.type() << " film model selected. "
67  << "Use the 'useReferenceValues' flag to employ reference "
68  << "pressure and temperature" << exit(FatalError);
69  }
70 
71  return refCast<const thermoSingleLayer>(filmModel_);
72 }
73 
74 
76 {
77  if (liquidPtr_ != nullptr)
78  {
79  return;
80  }
81 
82  dict.readEntry("liquid", name_);
83 
84  const SLGThermo* thermoPtr =
86 
87  if (thermoPtr)
88  {
89  // Retrieve from film thermo
90  ownLiquid_ = false;
91 
92  const SLGThermo& thermo = *thermoPtr;
93 
94  const label id = thermo.liquidId(name_);
95 
96  liquidPtr_ = &thermo.liquids().properties()[id];
97  }
98  else
99  {
100  // New liquid create
101  ownLiquid_ = true;
102 
103  liquidPtr_ =
104  liquidProperties::New(dict.optionalSubDict(name_ + "Coeffs")).ptr();
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  surfaceFilmRegionModel& film,
114  const dictionary& dict
115 )
116 :
117  filmThermoModel(typeName, film, dict),
118  name_("unknown_liquid"),
119  liquidPtr_(nullptr),
120  ownLiquid_(false),
121  useReferenceValues_(coeffDict_.get<bool>("useReferenceValues")),
122  pRef_(0.0),
123  TRef_(0.0)
124 {
126 
128  {
129  coeffDict_.readEntry("pRef", pRef_);
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 
138 {
139  if (ownLiquid_)
140  {
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
147 
149 {
150  return name_;
151 }
152 
153 
155 (
156  const scalar p,
157  const scalar T
158 ) const
159 {
160  return liquidPtr_->rho(p, T);
161 }
162 
163 
165 (
166  const scalar p,
167  const scalar T
168 ) const
169 {
170  return liquidPtr_->mu(p, T);
171 }
172 
173 
175 (
176  const scalar p,
177  const scalar T
178 ) const
179 {
180  return liquidPtr_->sigma(p, T);
181 }
182 
183 
185 (
186  const scalar p,
187  const scalar T
188 ) const
189 {
190  return liquidPtr_->Cp(p, T);
191 }
192 
193 
195 (
196  const scalar p,
197  const scalar T
198 ) const
199 {
200  return liquidPtr_->kappa(p, T);
201 }
202 
203 
204 scalar liquidFilmThermo::D
205 (
206  const scalar p,
207  const scalar T
208 ) const
209 {
210  return liquidPtr_->D(p, T);
211 }
212 
213 
215 (
216  const scalar p,
217  const scalar T
218 ) const
219 {
220  return liquidPtr_->hl(p, T);
221 }
222 
223 
225 (
226  const scalar p,
227  const scalar T
228 ) const
229 {
230  return liquidPtr_->pv(p, T);
231 }
232 
234 scalar liquidFilmThermo::W() const
235 {
236  return liquidPtr_->W();
237 }
238 
240 scalar liquidFilmThermo::Tb(const scalar p) const
241 {
242  return liquidPtr_->pvInvert(p);
243 }
244 
245 
247 {
249  (
250  IOobject::scopedName(type(), "rho"),
252  film().regionMesh(),
255  );
256  scalarField& rho = trho.ref().primitiveFieldRef();
257 
259  {
260  rho = this->rho(pRef_, TRef_);
261  }
262  else
263  {
264  const thermoSingleLayer& film = thermoFilm();
265 
266  const volScalarField& T = film.T();
267  const volScalarField& p = film.pPrimary();
268 
269  forAll(rho, celli)
270  {
271  rho[celli] = this->rho(p[celli], T[celli]);
272  }
273  }
275  trho.ref().correctBoundaryConditions();
276 
277  return trho;
278 }
279 
280 
282 {
283  auto tmu = volScalarField::New
284  (
285  IOobject::scopedName(type(), "mu"),
287  film().regionMesh(),
289  extrapolatedCalculatedFvPatchScalarField::typeName
290  );
291  scalarField& mu = tmu.ref().primitiveFieldRef();
292 
294  {
295  mu = this->mu(pRef_, TRef_);
296  }
297  else
298  {
299  const thermoSingleLayer& film = thermoFilm();
300 
301  const volScalarField& T = film.T();
302  const volScalarField& p = film.pPrimary();
303 
304  forAll(mu, celli)
305  {
306  mu[celli] = this->mu(p[celli], T[celli]);
307  }
308  }
310  tmu.ref().correctBoundaryConditions();
311 
312  return tmu;
313 }
314 
315 
317 {
318  auto tsigma = volScalarField::New
319  (
320  IOobject::scopedName(type(), "sigma"),
322  film().regionMesh(),
324  extrapolatedCalculatedFvPatchScalarField::typeName
325  );
326  scalarField& sigma = tsigma.ref().primitiveFieldRef();
327 
329  {
330  sigma = this->sigma(pRef_, TRef_);
331  }
332  else
333  {
334  const thermoSingleLayer& film = thermoFilm();
335 
336  const volScalarField& T = film.T();
337  const volScalarField& p = film.pPrimary();
338 
339  forAll(sigma, celli)
340  {
341  sigma[celli] = this->sigma(p[celli], T[celli]);
342  }
343  }
345  tsigma.ref().correctBoundaryConditions();
346 
347  return tsigma;
348 }
349 
350 
352 {
353  auto tCp = volScalarField::New
354  (
355  IOobject::scopedName(type(), "Cp"),
357  film().regionMesh(),
359  extrapolatedCalculatedFvPatchScalarField::typeName
360  );
361  scalarField& Cp = tCp.ref().primitiveFieldRef();
362 
364  {
365  Cp = this->Cp(pRef_, TRef_);
366  }
367  else
368  {
369  const thermoSingleLayer& film = thermoFilm();
370 
371  const volScalarField& T = film.T();
372  const volScalarField& p = film.pPrimary();
373 
374  forAll(Cp, celli)
375  {
376  Cp[celli] = this->Cp(p[celli], T[celli]);
377  }
378  }
380  tCp.ref().correctBoundaryConditions();
381 
382  return tCp;
383 }
384 
385 
387 {
388  auto tkappa = volScalarField::New
389  (
390  IOobject::scopedName(type(), "kappa"),
392  film().regionMesh(),
394  extrapolatedCalculatedFvPatchScalarField::typeName
395  );
396  scalarField& kappa = tkappa.ref().primitiveFieldRef();
397 
399  {
400  kappa = this->kappa(pRef_, TRef_);
401  }
402  else
403  {
404  const thermoSingleLayer& film = thermoFilm();
405 
406  const volScalarField& T = film.T();
407  const volScalarField& p = film.pPrimary();
408 
409  forAll(kappa, celli)
410  {
411  kappa[celli] = this->kappa(p[celli], T[celli]);
412  }
413  }
414 
415  tkappa.ref().correctBoundaryConditions();
416 
417  return tkappa;
418 }
419 
420 
421 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
422 
423 } // End namespace surfaceFilmModels
424 } // End namespace regionModels
425 } // End namespace Foam
426 
427 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
dictionary dict
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
virtual scalar rho(scalar p, scalar T) const =0
Liquid density [kg/m^3].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
tmp< volScalarField > trho
virtual scalar hl(scalar p, scalar T) const =0
Heat of vapourisation [J/kg].
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:92
liquidFilmThermo(const liquidFilmThermo &)=delete
No copy construct.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:104
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
virtual scalar kappa(scalar p, scalar T) const =0
Liquid thermal conductivity [W/(m K)].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual scalar sigma(scalar p, scalar T) const =0
Surface tension [N/m].
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
virtual scalar W() const
Return molecular weight [kg/kmol].
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
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:560
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
scalar W() const
Molecular weight [kg/kmol].
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
virtual scalar Cp(const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/(kg K)].
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
virtual scalar pv(scalar p, scalar T) const =0
Vapour pressure [Pa].
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:60
const dimensionSet dimPressure
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const dimensionSet dimPower
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
static autoPtr< liquidProperties > New(const word &name)
Return a pointer to a new liquidProperties created from name.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual scalar D(scalar p, scalar T) const =0
Vapour diffusivity [m2/s].
const dimensionSet dimEnergy
const dimensionSet dimDensity
virtual const word & name() const
Return the specie name.
Template functions to aid in the implementation of demand driven data.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
Definition: regionModel.H:239
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual scalar mu(scalar p, scalar T) const =0
Liquid viscosity [Pa s].
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
void deleteDemandDrivenData(DataPtr &dataPtr)
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Do not request registration (bool: false)
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
bool ownLiquid_
Flag to indicate that model owns the liquid object.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual const volScalarField & T() const =0
Return the film mean temperature [K].