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-2015 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"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 namespace surfaceFilmModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
44 
46 (
49  dictionary
50 );
51 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 {
57  if (!isA<thermoSingleLayer>(owner_))
58  {
60  << "Thermo model requires a " << thermoSingleLayer::typeName
61  << " film to supply the pressure and temperature, but "
62  << owner_.type() << " film model selected. "
63  << "Use the 'useReferenceValues' flag to employ reference "
64  << "pressure and temperature" << exit(FatalError);
65  }
66 
67  return refCast<const thermoSingleLayer>(owner_);
68 }
69 
70 
72 {
73  if (liquidPtr_ != NULL)
74  {
75  return;
76  }
77 
78  dict.lookup("liquid") >> name_;
79 
80  if (owner_.primaryMesh().foundObject<SLGThermo>("SLGThermo"))
81  {
82  // retrieve from film thermo
83  ownLiquid_ = false;
84 
85  const SLGThermo& thermo =
86  owner_.primaryMesh().lookupObject<SLGThermo>("SLGThermo");
87  label id = thermo.liquidId(name_);
88  liquidPtr_ = &thermo.liquids().properties()[id];
89  }
90  else
91  {
92  // new liquid create
93  ownLiquid_ = true;
94 
95  liquidPtr_ = new liquidProperties(dict.subDict(name_ + "Coeffs"));
96  }
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
103 (
104  surfaceFilmModel& owner,
105  const dictionary& dict
106 )
107 :
108  filmThermoModel(typeName, owner, dict),
109  name_("unknown_liquid"),
110  liquidPtr_(NULL),
111  ownLiquid_(false),
112  useReferenceValues_(readBool(coeffDict_.lookup("useReferenceValues"))),
113  pRef_(0.0),
114  TRef_(0.0)
115 {
116  initLiquid(coeffDict_);
117 
118  if (useReferenceValues_)
119  {
120  coeffDict_.lookup("pRef") >> pRef_;
121  coeffDict_.lookup("TRef") >> TRef_;
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 
129 {
130  if (ownLiquid_)
131  {
133  }
134 }
135 
136 
137 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
138 
140 {
141  return name_;
142 }
143 
144 
146 (
147  const scalar p,
148  const scalar T
149 ) const
150 {
151  return liquidPtr_->rho(p, T);
152 }
153 
154 
156 (
157  const scalar p,
158  const scalar T
159 ) const
160 {
161  return liquidPtr_->mu(p, T);
162 }
163 
164 
166 (
167  const scalar p,
168  const scalar T
169 ) const
170 {
171  return liquidPtr_->sigma(p, T);
172 }
173 
174 
176 (
177  const scalar p,
178  const scalar T
179 ) const
180 {
181  return liquidPtr_->Cp(p, T);
182 }
183 
184 
186 (
187  const scalar p,
188  const scalar T
189 ) const
190 {
191  return liquidPtr_->K(p, T);
192 }
193 
194 
195 scalar liquidFilmThermo::D
196 (
197  const scalar p,
198  const scalar T
199 ) const
200 {
201  return liquidPtr_->D(p, T);
202 }
203 
204 
206 (
207  const scalar p,
208  const scalar T
209 ) const
210 {
211  return liquidPtr_->hl(p, T);
212 }
213 
214 
216 (
217  const scalar p,
218  const scalar T
219 ) const
220 {
221  return liquidPtr_->pv(p, T);
222 }
223 
224 
225 scalar liquidFilmThermo::W() const
226 {
227  return liquidPtr_->W();
228 }
229 
230 
231 scalar liquidFilmThermo::Tb(const scalar p) const
232 {
233  return liquidPtr_->pvInvert(p);
234 }
235 
236 
238 {
240  (
241  new volScalarField
242  (
243  IOobject
244  (
245  type() + ":rho",
246  owner().time().timeName(),
247  owner().regionMesh(),
250  ),
251  owner().regionMesh(),
252  dimensionedScalar("0", dimDensity, 0.0),
253  zeroGradientFvPatchScalarField::typeName
254  )
255  );
256 
257  scalarField& rho = trho().internalField();
258 
260  {
261  rho = this->rho(pRef_, TRef_);
262  }
263  else
264  {
265  const thermoSingleLayer& film = thermoFilm();
266 
267  const volScalarField& T = film.T();
268  const volScalarField& p = film.pPrimary();
269 
270  forAll(rho, cellI)
271  {
272  rho[cellI] = this->rho(p[cellI], T[cellI]);
273  }
274  }
275 
276  trho().correctBoundaryConditions();
277 
278  return trho;
279 }
280 
281 
283 {
285  (
286  new volScalarField
287  (
288  IOobject
289  (
290  type() + ":mu",
291  owner().time().timeName(),
292  owner().regionMesh(),
295  ),
296  owner().regionMesh(),
298  zeroGradientFvPatchScalarField::typeName
299  )
300  );
301 
302  scalarField& mu = tmu().internalField();
303 
305  {
306  mu = this->mu(pRef_, TRef_);
307  }
308  else
309  {
310  const thermoSingleLayer& film = thermoFilm();
311 
312  const volScalarField& T = film.T();
313  const volScalarField& p = film.pPrimary();
314 
315  forAll(mu, cellI)
316  {
317  mu[cellI] = this->mu(p[cellI], T[cellI]);
318  }
319  }
320 
321  tmu().correctBoundaryConditions();
322 
323  return tmu;
324 }
325 
326 
328 {
329  tmp<volScalarField> tsigma
330  (
331  new volScalarField
332  (
333  IOobject
334  (
335  type() + ":sigma",
336  owner().time().timeName(),
337  owner().regionMesh(),
340  ),
341  owner().regionMesh(),
342  dimensionedScalar("0", dimMass/sqr(dimTime), 0.0),
343  zeroGradientFvPatchScalarField::typeName
344  )
345  );
346 
347  scalarField& sigma = tsigma().internalField();
348 
350  {
351  sigma = this->sigma(pRef_, TRef_);
352  }
353  else
354  {
355  const thermoSingleLayer& film = thermoFilm();
356 
357  const volScalarField& T = film.T();
358  const volScalarField& p = film.pPrimary();
359 
360  forAll(sigma, cellI)
361  {
362  sigma[cellI] = this->sigma(p[cellI], T[cellI]);
363  }
364  }
365 
366  tsigma().correctBoundaryConditions();
367 
368  return tsigma;
369 }
370 
371 
373 {
375  (
376  new volScalarField
377  (
378  IOobject
379  (
380  type() + ":Cp",
381  owner().time().timeName(),
382  owner().regionMesh(),
385  ),
386  owner().regionMesh(),
388  zeroGradientFvPatchScalarField::typeName
389  )
390  );
391 
392  scalarField& Cp = tCp().internalField();
393 
395  {
396  Cp = this->Cp(pRef_, TRef_);
397  }
398  else
399  {
400  const thermoSingleLayer& film = thermoFilm();
401 
402  const volScalarField& T = film.T();
403  const volScalarField& p = film.pPrimary();
404 
405  forAll(Cp, cellI)
406  {
407  Cp[cellI] = this->Cp(p[cellI], T[cellI]);
408  }
409  }
410 
411  tCp().correctBoundaryConditions();
412 
413  return tCp;
414 }
415 
416 
418 {
419  tmp<volScalarField> tkappa
420  (
421  new volScalarField
422  (
423  IOobject
424  (
425  type() + ":kappa",
426  owner().time().timeName(),
427  owner().regionMesh(),
430  ),
431  owner().regionMesh(),
433  zeroGradientFvPatchScalarField::typeName
434  )
435  );
436 
437  scalarField& kappa = tkappa().internalField();
438 
440  {
441  kappa = this->kappa(pRef_, TRef_);
442  }
443  else
444  {
445  const thermoSingleLayer& film = thermoFilm();
446 
447  const volScalarField& T = film.T();
448  const volScalarField& p = film.pPrimary();
449 
450  forAll(kappa, cellI)
451  {
452  kappa[cellI] = this->kappa(p[cellI], T[cellI]);
453  }
454  }
455 
456  tkappa().correctBoundaryConditions();
457 
458  return tkappa;
459 }
460 
461 
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
463 
464 } // End namespace surfaceFilmModels
465 } // End namespace regionModels
466 } // End namespace Foam
467 
468 // ************************************************************************* //
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::thermoFilm
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
Definition: liquidFilmThermo.C:55
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::pRef_
scalar pRef_
Reference pressure [pa].
Definition: liquidFilmThermo.H:77
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
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:139
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::kappa
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: liquidFilmThermo.C:417
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:196
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::sigma
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
Definition: liquidFilmThermo.C:327
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
Definition: filmSubModelBaseI.H:37
SLGThermo.H
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::rho
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
Definition: liquidFilmThermo.C:237
Foam::liquidProperties::W
scalar W() const
Molecular weight [kg/kmol].
Definition: liquidPropertiesI.H:26
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
thermoSingleLayer.H
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::liquidFilmThermo
liquidFilmThermo(const liquidFilmThermo &)
Disallow default bitwise copy construct.
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::ownLiquid_
bool ownLiquid_
Flag to indicate that model owns the liquid object.
Definition: liquidFilmThermo.H:71
Foam::regionModels::regionModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
liquidFilmThermo
Liquid thermo model.
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
liquidFilmThermo.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::liquidProperties
The thermophysical properties of a liquidProperties.
Definition: liquidProperties.H:53
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Definition: thermoSingleLayer.H:63
Foam::label
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
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::initLiquid
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
Definition: liquidFilmThermo.C:71
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::Cp
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: liquidFilmThermo.C:372
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::subModelBase::dict
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::mu
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
Definition: liquidFilmThermo.C:282
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::liquidPtr_
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
Definition: liquidFilmThermo.H:68
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary
const volScalarField & pPrimary() const
Pressure / [Pa].
Definition: kinematicSingleLayerI.H:155
Foam::liquidProperties::pvInvert
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
Definition: liquidProperties.C:314
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: thermoSingleLayer.C:707
Foam::dimPower
const dimensionSet dimPower
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
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:206
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::Tb
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
Definition: liquidFilmThermo.C:231
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::useReferenceValues_
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
Definition: liquidFilmThermo.H:74
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::regionModels::surfaceFilmModels::filmThermoModel
Definition: filmThermoModel.H:54
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
filmThermoModel
Base class for film thermo models.
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::TRef_
scalar TRef_
Reference temperature [K].
Definition: liquidFilmThermo.H:80
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner_
surfaceFilmModel & owner_
Reference to the owner surface film model.
Definition: filmSubModelBase.H:63
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::W
virtual scalar W() const
Return molecular weight [kg/kmol].
Definition: liquidFilmThermo.C:225
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::name_
word name_
Liquid name.
Definition: liquidFilmThermo.H:65
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::~liquidFilmThermo
virtual ~liquidFilmThermo()
Destructor.
Definition: liquidFilmThermo.C:128
Foam::regionModels::surfaceFilmModels::liquidFilmThermo::pv
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
Definition: liquidFilmThermo.C:216