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-2016 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>(owner_))
59  {
61  << "Thermo model requires a " << thermoSingleLayer::typeName
62  << " film to supply the pressure and temperature, but "
63  << owner_.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>(owner_);
69 }
70 
71 
73 {
74  if (liquidPtr_ != NULL)
75  {
76  return;
77  }
78 
79  dict.lookup("liquid") >> name_;
80 
81  if (owner_.primaryMesh().foundObject<SLGThermo>("SLGThermo"))
82  {
83  // retrieve from film thermo
84  ownLiquid_ = false;
85 
86  const SLGThermo& thermo =
87  owner_.primaryMesh().lookupObject<SLGThermo>("SLGThermo");
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_ = new liquidProperties(dict.subDict(name_ + "Coeffs"));
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 (
105  surfaceFilmModel& owner,
106  const dictionary& dict
107 )
108 :
109  filmThermoModel(typeName, owner, dict),
110  name_("unknown_liquid"),
111  liquidPtr_(NULL),
112  ownLiquid_(false),
113  useReferenceValues_(readBool(coeffDict_.lookup("useReferenceValues"))),
114  pRef_(0.0),
115  TRef_(0.0)
116 {
118 
120  {
121  coeffDict_.lookup("pRef") >> pRef_;
122  coeffDict_.lookup("TRef") >> TRef_;
123  }
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
128 
130 {
131  if (ownLiquid_)
132  {
134  }
135 }
136 
137 
138 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
139 
141 {
142  return name_;
143 }
144 
145 
147 (
148  const scalar p,
149  const scalar T
150 ) const
151 {
152  return liquidPtr_->rho(p, T);
153 }
154 
155 
157 (
158  const scalar p,
159  const scalar T
160 ) const
161 {
162  return liquidPtr_->mu(p, T);
163 }
164 
165 
167 (
168  const scalar p,
169  const scalar T
170 ) const
171 {
172  return liquidPtr_->sigma(p, T);
173 }
174 
175 
177 (
178  const scalar p,
179  const scalar T
180 ) const
181 {
182  return liquidPtr_->Cp(p, T);
183 }
184 
185 
187 (
188  const scalar p,
189  const scalar T
190 ) const
191 {
192  return liquidPtr_->K(p, T);
193 }
194 
195 
196 scalar liquidFilmThermo::D
197 (
198  const scalar p,
199  const scalar T
200 ) const
201 {
202  return liquidPtr_->D(p, T);
203 }
204 
205 
207 (
208  const scalar p,
209  const scalar T
210 ) const
211 {
212  return liquidPtr_->hl(p, T);
213 }
214 
215 
217 (
218  const scalar p,
219  const scalar T
220 ) const
221 {
222  return liquidPtr_->pv(p, T);
223 }
224 
225 
226 scalar liquidFilmThermo::W() const
227 {
228  return liquidPtr_->W();
229 }
230 
231 
232 scalar liquidFilmThermo::Tb(const scalar p) const
233 {
234  return liquidPtr_->pvInvert(p);
235 }
236 
237 
239 {
241  (
242  new volScalarField
243  (
244  IOobject
245  (
246  type() + ":rho",
247  owner().time().timeName(),
248  owner().regionMesh(),
251  ),
252  owner().regionMesh(),
253  dimensionedScalar("0", dimDensity, 0.0),
254  extrapolatedCalculatedFvPatchScalarField::typeName
255  )
256  );
257 
258  scalarField& rho = trho.ref().primitiveFieldRef();
259 
261  {
262  forAll(rho, celli)
263  {
264  rho[celli] = this->rho(pRef_, TRef_);
265  }
266  }
267  else
268  {
269  const thermoSingleLayer& film = thermoFilm();
270 
271  const volScalarField& T = film.T();
272  const volScalarField& p = film.pPrimary();
273 
274  forAll(rho, celli)
275  {
276  rho[celli] = this->rho(p[celli], T[celli]);
277  }
278  }
279 
280  trho.ref().correctBoundaryConditions();
281 
282  return trho;
283 }
284 
285 
287 {
289  (
290  new volScalarField
291  (
292  IOobject
293  (
294  type() + ":mu",
295  owner().time().timeName(),
296  owner().regionMesh(),
299  ),
300  owner().regionMesh(),
302  extrapolatedCalculatedFvPatchScalarField::typeName
303  )
304  );
305 
306  scalarField& mu = tmu.ref().primitiveFieldRef();
307 
309  {
310  forAll(mu, celli)
311  {
312  mu[celli] = this->mu(pRef_, TRef_);
313  }
314  }
315  else
316  {
317  const thermoSingleLayer& film = thermoFilm();
318 
319  const volScalarField& T = film.T();
320  const volScalarField& p = film.pPrimary();
321 
322  forAll(mu, celli)
323  {
324  mu[celli] = this->mu(p[celli], T[celli]);
325  }
326  }
327 
328  tmu.ref().correctBoundaryConditions();
329 
330  return tmu;
331 }
332 
333 
335 {
336  tmp<volScalarField> tsigma
337  (
338  new volScalarField
339  (
340  IOobject
341  (
342  type() + ":sigma",
343  owner().time().timeName(),
344  owner().regionMesh(),
347  ),
348  owner().regionMesh(),
349  dimensionedScalar("0", dimMass/sqr(dimTime), 0.0),
350  extrapolatedCalculatedFvPatchScalarField::typeName
351  )
352  );
353 
354  scalarField& sigma = tsigma.ref().primitiveFieldRef();
355 
357  {
358  forAll(sigma, celli)
359  {
360  sigma[celli] = this->sigma(pRef_, TRef_);
361  }
362  }
363  else
364  {
365  const thermoSingleLayer& film = thermoFilm();
366 
367  const volScalarField& T = film.T();
368  const volScalarField& p = film.pPrimary();
369 
370  forAll(sigma, celli)
371  {
372  sigma[celli] = this->sigma(p[celli], T[celli]);
373  }
374  }
375 
376  tsigma.ref().correctBoundaryConditions();
377 
378  return tsigma;
379 }
380 
381 
383 {
385  (
386  new volScalarField
387  (
388  IOobject
389  (
390  type() + ":Cp",
391  owner().time().timeName(),
392  owner().regionMesh(),
395  ),
396  owner().regionMesh(),
398  extrapolatedCalculatedFvPatchScalarField::typeName
399  )
400  );
401 
402  scalarField& Cp = tCp.ref().primitiveFieldRef();
403 
405  {
406  forAll(Cp, celli)
407  {
408  Cp[celli] = this->Cp(pRef_, TRef_);
409  }
410  }
411  else
412  {
413  const thermoSingleLayer& film = thermoFilm();
414 
415  const volScalarField& T = film.T();
416  const volScalarField& p = film.pPrimary();
417 
418  forAll(Cp, celli)
419  {
420  Cp[celli] = this->Cp(p[celli], T[celli]);
421  }
422  }
423 
424  tCp.ref().correctBoundaryConditions();
425 
426  return tCp;
427 }
428 
429 
431 {
432  tmp<volScalarField> tkappa
433  (
434  new volScalarField
435  (
436  IOobject
437  (
438  type() + ":kappa",
439  owner().time().timeName(),
440  owner().regionMesh(),
443  ),
444  owner().regionMesh(),
446  extrapolatedCalculatedFvPatchScalarField::typeName
447  )
448  );
449 
450  scalarField& kappa = tkappa.ref().primitiveFieldRef();
451 
453  {
454  forAll(kappa, celli)
455  {
456  kappa[celli] = this->kappa(pRef_, TRef_);
457  }
458  }
459  else
460  {
461  const thermoSingleLayer& film = thermoFilm();
462 
463  const volScalarField& T = film.T();
464  const volScalarField& p = film.pPrimary();
465 
466  forAll(kappa, celli)
467  {
468  kappa[celli] = this->kappa(p[celli], T[celli]);
469  }
470  }
471 
472  tkappa.ref().correctBoundaryConditions();
473 
474  return tkappa;
475 }
476 
477 
478 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
479 
480 } // End namespace surfaceFilmModels
481 } // End namespace regionModels
482 } // End namespace Foam
483 
484 // ************************************************************************* //
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
const PtrList< liquidProperties > & properties() const
Return the liquid properties.
label liquidId(const word &cmptName, bool allowNotFound=false) const
Index of liquid component.
Definition: SLGThermo.C:174
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual scalar W() const
Return molecular weight [kg/kmol].
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool foundObject(const word &name) const
Is the named Type found?
surfaceFilmModel & owner_
Reference to the owner surface film model.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
virtual scalar hl(scalar p, scalar T) const
Heat of vapourisation [J/kg].
virtual scalar mu(scalar p, scalar T) const
Liquid viscosity [Pa s].
tmp< volScalarField > trho
virtual const word & name() const
Return the specie name.
bool readBool(Istream &)
Definition: boolIO.C:60
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 > sigma() const
Return surface tension [kg/s2].
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
A class for handling words, derived from string.
Definition: word.H:59
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual scalar sigma(scalar p, scalar T) const
Surface tension [N/m].
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
The thermophysical properties of a liquidProperties.
word timeName
Definition: getTimeIndex.H:3
virtual scalar Cp(scalar p, scalar T) const
Liquid heat capacity [J/(kg K)].
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
const liquidMixtureProperties & liquids() const
Return reference to the global (additional) liquids.
Definition: SLGThermo.C:121
const dimensionSet dimPressure
virtual scalar rho(scalar p, scalar T) const
Liquid rho [kg/m^3].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const dimensionSet dimPower
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual scalar pv(scalar p, scalar T) const
Vapour pressure [Pa].
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
virtual scalar K(scalar p, scalar T) const
Liquid thermal conductivity [W/(m K)].
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
const dimensionSet dimEnergy
const dimensionSet dimDensity
Template functions to aid in the implementation of demand driven data.
scalar W() const
Molecular weight [kg/kmol].
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
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 surfaceFilmModel & owner() const
Return const access to the owner surface film model.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
virtual scalar D(scalar p, scalar T) const
Vapour diffussivity [m2/s].
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
void deleteDemandDrivenData(DataPtr &dataPtr)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
defineTypeNameAndDebug(kinematicSingleLayer, 0)
const volScalarField & pPrimary() const
Pressure / [Pa].
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Namespace for OpenFOAM.
bool ownLiquid_
Flag to indicate that model owns the liquid object.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451