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 | Copyright (C) 2013-2017 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 (
106  surfaceFilmModel& film,
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  (
243  new volScalarField
244  (
245  IOobject
246  (
247  type() + ":rho",
248  film().time().timeName(),
249  film().regionMesh(),
252  ),
253  film().regionMesh(),
254  dimensionedScalar("0", dimDensity, 0.0),
255  extrapolatedCalculatedFvPatchScalarField::typeName
256  )
257  );
258 
259  scalarField& rho = trho.ref().primitiveFieldRef();
260 
262  {
263  forAll(rho, celli)
264  {
265  rho[celli] = this->rho(pRef_, TRef_);
266  }
267  }
268  else
269  {
270  const thermoSingleLayer& film = thermoFilm();
271 
272  const volScalarField& T = film.T();
273  const volScalarField& p = film.pPrimary();
274 
275  forAll(rho, celli)
276  {
277  rho[celli] = this->rho(p[celli], T[celli]);
278  }
279  }
280 
281  trho.ref().correctBoundaryConditions();
282 
283  return trho;
284 }
285 
286 
288 {
290  (
291  new volScalarField
292  (
293  IOobject
294  (
295  type() + ":mu",
296  film().time().timeName(),
297  film().regionMesh(),
300  ),
301  film().regionMesh(),
303  extrapolatedCalculatedFvPatchScalarField::typeName
304  )
305  );
306 
307  scalarField& mu = tmu.ref().primitiveFieldRef();
308 
310  {
311  forAll(mu, celli)
312  {
313  mu[celli] = this->mu(pRef_, TRef_);
314  }
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(),
350  dimensionedScalar("0", dimMass/sqr(dimTime), 0.0),
351  extrapolatedCalculatedFvPatchScalarField::typeName
352  )
353  );
354 
355  scalarField& sigma = tsigma.ref().primitiveFieldRef();
356 
358  {
359  forAll(sigma, celli)
360  {
361  sigma[celli] = this->sigma(pRef_, TRef_);
362  }
363  }
364  else
365  {
366  const thermoSingleLayer& film = thermoFilm();
367 
368  const volScalarField& T = film.T();
369  const volScalarField& p = film.pPrimary();
370 
371  forAll(sigma, celli)
372  {
373  sigma[celli] = this->sigma(p[celli], T[celli]);
374  }
375  }
376 
377  tsigma.ref().correctBoundaryConditions();
378 
379  return tsigma;
380 }
381 
382 
384 {
386  (
387  new volScalarField
388  (
389  IOobject
390  (
391  type() + ":Cp",
392  film().time().timeName(),
393  film().regionMesh(),
396  ),
397  film().regionMesh(),
399  extrapolatedCalculatedFvPatchScalarField::typeName
400  )
401  );
402 
403  scalarField& Cp = tCp.ref().primitiveFieldRef();
404 
406  {
407  forAll(Cp, celli)
408  {
409  Cp[celli] = this->Cp(pRef_, TRef_);
410  }
411  }
412  else
413  {
414  const thermoSingleLayer& film = thermoFilm();
415 
416  const volScalarField& T = film.T();
417  const volScalarField& p = film.pPrimary();
418 
419  forAll(Cp, celli)
420  {
421  Cp[celli] = this->Cp(p[celli], T[celli]);
422  }
423  }
424 
425  tCp.ref().correctBoundaryConditions();
426 
427  return tCp;
428 }
429 
430 
432 {
433  tmp<volScalarField> tkappa
434  (
435  new volScalarField
436  (
437  IOobject
438  (
439  type() + ":kappa",
440  film().time().timeName(),
441  film().regionMesh(),
444  ),
445  film().regionMesh(),
447  extrapolatedCalculatedFvPatchScalarField::typeName
448  )
449  );
450 
451  scalarField& kappa = tkappa.ref().primitiveFieldRef();
452 
454  {
455  forAll(kappa, celli)
456  {
457  kappa[celli] = this->kappa(pRef_, TRef_);
458  }
459  }
460  else
461  {
462  const thermoSingleLayer& film = thermoFilm();
463 
464  const volScalarField& T = film.T();
465  const volScalarField& p = film.pPrimary();
466 
467  forAll(kappa, celli)
468  {
469  kappa[celli] = this->kappa(p[celli], T[celli]);
470  }
471  }
472 
473  tkappa.ref().correctBoundaryConditions();
474 
475  return tkappa;
476 }
477 
478 
479 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480 
481 } // End namespace surfaceFilmModels
482 } // End namespace regionModels
483 } // End namespace Foam
484 
485 // ************************************************************************* //
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
const surfaceFilmModel & film() const
Return const access to the film surface film model.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:137
#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
surfaceFilmModel & filmModel_
Reference to the film surface film model.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:79
liquidFilmThermo(const liquidFilmThermo &)
Disallow default bitwise copy construct.
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)].
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.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
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:759
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/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].
word timeName
Definition: getTimeIndex.H:3
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/s2].
const volScalarField & pPrimary() const
Pressure [Pa].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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].
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
virtual const volScalarField & T() const
Return the film mean temperature [K].
bool ownLiquid_
Flag to indicate that model owns the liquid object.