liquidFilmThermo.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) 2013-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 "liquidFilmThermo.H"
27 #include "demandDrivenData.H"
28 #include "thermoSingleLayer.H"
29 #include "SLGThermo.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(liquidFilmThermo, 0);
45 
47 (
48  filmThermoModel,
49  liquidFilmThermo,
50  dictionary
51 );
52 
53 
54 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55 
57 {
58  if (!isA<thermoSingleLayer>(filmModel_))
59  {
61  << "Thermo model requires a " << thermoSingleLayer::typeName
62  << " film to supply the pressure and temperature, but "
63  << filmModel_.type() << " film model selected. "
64  << "Use the 'useReferenceValues' flag to employ reference "
65  << "pressure and temperature" << exit(FatalError);
66  }
67 
68  return refCast<const thermoSingleLayer>(filmModel_);
69 }
70 
71 
73 {
74  if (liquidPtr_ != nullptr)
75  {
76  return;
77  }
78 
79  dict.lookup("liquid") >> name_;
80 
81  if (filmModel_.primaryMesh().foundObject<SLGThermo>("SLGThermo"))
82  {
83  // retrieve from film thermo
84  ownLiquid_ = false;
85 
86  const SLGThermo& thermo =
88  label id = thermo.liquidId(name_);
89  liquidPtr_ = &thermo.liquids().properties()[id];
90  }
91  else
92  {
93  // new liquid create
94  ownLiquid_ = true;
95 
96  liquidPtr_ =
97  liquidProperties::New(dict.optionalSubDict(name_ + "Coeffs")).ptr();
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
105 (
107  const dictionary& dict
108 )
109 :
110  filmThermoModel(typeName, film, dict),
111  name_("unknown_liquid"),
112  liquidPtr_(nullptr),
113  ownLiquid_(false),
114  useReferenceValues_(readBool(coeffDict_.lookup("useReferenceValues"))),
115  pRef_(0.0),
116  TRef_(0.0)
117 {
119 
121  {
122  coeffDict_.lookup("pRef") >> pRef_;
123  coeffDict_.lookup("TRef") >> TRef_;
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
131 {
132  if (ownLiquid_)
133  {
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140 
142 {
143  return name_;
144 }
145 
146 
148 (
149  const scalar p,
150  const scalar T
151 ) const
152 {
153  return liquidPtr_->rho(p, T);
154 }
155 
156 
158 (
159  const scalar p,
160  const scalar T
161 ) const
162 {
163  return liquidPtr_->mu(p, T);
164 }
165 
166 
168 (
169  const scalar p,
170  const scalar T
171 ) const
172 {
173  return liquidPtr_->sigma(p, T);
174 }
175 
176 
178 (
179  const scalar p,
180  const scalar T
181 ) const
182 {
183  return liquidPtr_->Cp(p, T);
184 }
185 
186 
188 (
189  const scalar p,
190  const scalar T
191 ) const
192 {
193  return liquidPtr_->kappa(p, T);
194 }
195 
196 
197 scalar liquidFilmThermo::D
198 (
199  const scalar p,
200  const scalar T
201 ) const
202 {
203  return liquidPtr_->D(p, T);
204 }
205 
206 
208 (
209  const scalar p,
210  const scalar T
211 ) const
212 {
213  return liquidPtr_->hl(p, T);
214 }
215 
216 
218 (
219  const scalar p,
220  const scalar T
221 ) const
222 {
223  return liquidPtr_->pv(p, T);
224 }
225 
226 
227 scalar liquidFilmThermo::W() const
228 {
229  return liquidPtr_->W();
230 }
231 
232 
233 scalar liquidFilmThermo::Tb(const scalar p) const
234 {
235  return liquidPtr_->pvInvert(p);
236 }
237 
238 
240 {
242  (
244  (
245  type() + ":rho",
246  film().regionMesh(),
248  extrapolatedCalculatedFvPatchScalarField::typeName
249  )
250  );
251 
252  scalarField& rho = trho.ref().primitiveFieldRef();
253 
255  {
256  forAll(rho, celli)
257  {
258  rho[celli] = this->rho(pRef_, TRef_);
259  }
260  }
261  else
262  {
263  const thermoSingleLayer& film = thermoFilm();
264 
265  const volScalarField& T = film.T();
266  const volScalarField& p = film.pPrimary();
267 
268  forAll(rho, celli)
269  {
270  rho[celli] = this->rho(p[celli], T[celli]);
271  }
272  }
273 
274  trho.ref().correctBoundaryConditions();
275 
276  return trho;
277 }
278 
279 
281 {
283  (
285  (
286  type() + ":mu",
287  film().regionMesh(),
289  extrapolatedCalculatedFvPatchScalarField::typeName
290  )
291  );
292 
293  scalarField& mu = tmu.ref().primitiveFieldRef();
294 
296  {
297  forAll(mu, celli)
298  {
299  mu[celli] = this->mu(pRef_, TRef_);
300  }
301  }
302  else
303  {
304  const thermoSingleLayer& film = thermoFilm();
305 
306  const volScalarField& T = film.T();
307  const volScalarField& p = film.pPrimary();
308 
309  forAll(mu, celli)
310  {
311  mu[celli] = this->mu(p[celli], T[celli]);
312  }
313  }
314 
315  tmu.ref().correctBoundaryConditions();
316 
317  return tmu;
318 }
319 
320 
322 {
323  tmp<volScalarField> tsigma
324  (
326  (
327  type() + ":sigma",
328  film().regionMesh(),
330  extrapolatedCalculatedFvPatchScalarField::typeName
331  )
332  );
333 
334  scalarField& sigma = tsigma.ref().primitiveFieldRef();
335 
337  {
338  forAll(sigma, celli)
339  {
340  sigma[celli] = this->sigma(pRef_, TRef_);
341  }
342  }
343  else
344  {
345  const thermoSingleLayer& film = thermoFilm();
346 
347  const volScalarField& T = film.T();
348  const volScalarField& p = film.pPrimary();
349 
350  forAll(sigma, celli)
351  {
352  sigma[celli] = this->sigma(p[celli], T[celli]);
353  }
354  }
355 
356  tsigma.ref().correctBoundaryConditions();
357 
358  return tsigma;
359 }
360 
361 
363 {
365  (
367  (
368  type() + ":Cp",
369  film().regionMesh(),
371  extrapolatedCalculatedFvPatchScalarField::typeName
372  )
373  );
374 
375  scalarField& Cp = tCp.ref().primitiveFieldRef();
376 
378  {
379  forAll(Cp, celli)
380  {
381  Cp[celli] = this->Cp(pRef_, TRef_);
382  }
383  }
384  else
385  {
386  const thermoSingleLayer& film = thermoFilm();
387 
388  const volScalarField& T = film.T();
389  const volScalarField& p = film.pPrimary();
390 
391  forAll(Cp, celli)
392  {
393  Cp[celli] = this->Cp(p[celli], T[celli]);
394  }
395  }
396 
397  tCp.ref().correctBoundaryConditions();
398 
399  return tCp;
400 }
401 
402 
404 {
405  tmp<volScalarField> tkappa
406  (
408  (
409  type() + ":kappa",
410  film().regionMesh(),
412  extrapolatedCalculatedFvPatchScalarField::typeName
413  )
414  );
415 
416  scalarField& kappa = tkappa.ref().primitiveFieldRef();
417 
419  {
420  forAll(kappa, celli)
421  {
422  kappa[celli] = this->kappa(pRef_, TRef_);
423  }
424  }
425  else
426  {
427  const thermoSingleLayer& film = thermoFilm();
428 
429  const volScalarField& T = film.T();
430  const volScalarField& p = film.pPrimary();
431 
432  forAll(kappa, celli)
433  {
434  kappa[celli] = this->kappa(p[celli], T[celli]);
435  }
436  }
437 
438  tkappa.ref().correctBoundaryConditions();
439 
440  return tkappa;
441 }
442 
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 } // End namespace surfaceFilmModels
447 } // End namespace regionModels
448 } // End namespace Foam
449 
450 // ************************************************************************* //
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
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
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:124
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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
tmp< volScalarField > trho
bool foundObject(const word &name) const
Is the named Type found?
virtual scalar hl(scalar p, scalar T) const =0
Heat of vapourisation [J/kg].
bool readBool(Istream &)
Definition: boolIO.C:60
rhoReactionThermo & thermo
Definition: createFields.H:28
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:74
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:766
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
scalar W() const
Molecular weight [kg/kmol].
const liquidMixtureProperties & liquids() const
Return reference to the global (additional) liquids.
Definition: SLGThermo.C:121
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
virtual tmp< volScalarField > rho() const
Return density [kg/m^3].
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].
filmThermoModel(surfaceFilmRegionModel &film)
Construct null.
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:62
const dimensionSet dimPressure
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s^2].
const volScalarField & pPrimary() const
Pressure [Pa].
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.
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual scalar D(scalar p, scalar T) const =0
Vapour diffussivity [m2/s].
const dimensionSet dimEnergy
const dimensionSet dimDensity
label liquidId(const word &cmptName, bool allowNotFound=false) const
Index of liquid component.
Definition: SLGThermo.C:174
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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
A class for managing temporary objects.
Definition: PtrList.H:53
const PtrList< liquidProperties > & properties() const
Return the liquid properties.
virtual scalar mu(scalar p, scalar T) const =0
Liquid viscosity [Pa s].
void deleteDemandDrivenData(DataPtr &dataPtr)
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
liquidFilmThermo(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model and dictionary.
virtual const volScalarField & T() const
Return the film mean temperature [K].
bool ownLiquid_
Flag to indicate that model owns the liquid object.