constantFilmThermo.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 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 "constantFilmThermo.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace regionModels
34 {
35 namespace surfaceFilmModels
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
41 
43 (
46  dictionary
47 );
48 
49 
51 {
53  {
54  td.set_ = true;
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
62 (
63  surfaceFilmModel& owner,
64  const dictionary& dict
65 )
66 :
67  filmThermoModel(typeName, owner, dict),
68  name_(coeffDict_.lookup("specieName")),
69  rho0_("rho0"),
70  mu0_("mu0"),
71  sigma0_("sigma0"),
72  Cp0_("Cp0"),
73  kappa0_("kappa0"),
74  hl0_("hl0"),
75  pv0_("pv0"),
76  W0_("W0"),
77  Tb0_("Tb0")
78 {
79  init(rho0_);
80  init(mu0_);
81  init(sigma0_);
82  init(Cp0_);
83  init(kappa0_);
84  init(hl0_);
85  init(pv0_);
86  init(W0_);
87  init(Tb0_);
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
98 
100 {
101  return name_;
102 }
103 
104 
106 (
107  const scalar p,
108  const scalar T
109 ) const
110 {
111  if (!rho0_.set_)
112  {
113  coeffDict_.lookup(rho0_.name_) >> rho0_.value_;
114  rho0_.set_ = true;
115  }
116 
117  return rho0_.value_;
118 }
119 
120 
122 (
123  const scalar p,
124  const scalar T
125 ) const
126 {
127  if (!mu0_.set_)
128  {
129  coeffDict_.lookup(mu0_.name_) >> mu0_.value_;
130  mu0_.set_ = true;
131  }
132 
133  return mu0_.value_;
134 }
135 
136 
138 (
139  const scalar p,
140  const scalar T
141 ) const
142 {
143  if (!sigma0_.set_)
144  {
145  coeffDict_.lookup(sigma0_.name_) >> sigma0_.value_;
146  sigma0_.set_ = true;
147  }
148 
149  return sigma0_.value_;
150 }
151 
152 
154 (
155  const scalar p,
156  const scalar T
157 ) const
158 {
159  if (!Cp0_.set_)
160  {
161  coeffDict_.lookup(Cp0_.name_) >> Cp0_.value_;
162  Cp0_.set_ = true;
163  }
164 
165  return Cp0_.value_;
166 }
167 
168 
170 (
171  const scalar p,
172  const scalar T
173 ) const
174 {
175  if (!kappa0_.set_)
176  {
177  coeffDict_.lookup(kappa0_.name_) >> kappa0_.value_;
178  kappa0_.set_ = true;
179  }
180 
181  return kappa0_.value_;
182 }
183 
184 
186 (
187  const scalar p,
188  const scalar T
189 ) const
190 {
191  if (!D0_.set_)
192  {
193  coeffDict_.lookup(D0_.name_) >> D0_.value_;
194  D0_.set_ = true;
195  }
196 
197  return D0_.value_;
198 }
199 
200 
202 (
203  const scalar p,
204  const scalar T
205 ) const
206 {
207  if (!hl0_.set_)
208  {
209  coeffDict_.lookup(hl0_.name_) >> hl0_.value_;
210  hl0_.set_ = true;
211  }
212 
213  return hl0_.value_;
214 }
215 
216 
218 (
219  const scalar p,
220  const scalar T
221 ) const
222 {
223  if (!pv0_.set_)
224  {
225  coeffDict_.lookup(pv0_.name_) >> pv0_.value_;
226  pv0_.set_ = true;
227  }
228 
229  return pv0_.value_;
230 }
231 
232 
233 scalar constantFilmThermo::W() const
234 {
235  if (!W0_.set_)
236  {
238  W0_.set_ = true;
239  }
240 
241  return W0_.value_;
242 }
243 
244 
245 scalar constantFilmThermo::Tb(const scalar p) const
246 {
247  if (!Tb0_.set_)
248  {
250  Tb0_.set_ = true;
251  }
252 
253  return Tb0_.value_;
254 }
255 
256 
258 {
260  (
261  new volScalarField
262  (
263  IOobject
264  (
265  type() + ':' + rho0_.name_,
266  owner().time().timeName(),
267  owner().regionMesh(),
270  ),
271  owner().regionMesh(),
272  dimensionedScalar("0", dimDensity, 0.0),
273  zeroGradientFvPatchScalarField::typeName
274  )
275  );
276 
277  trho().internalField() = this->rho(0, 0);
278  trho().correctBoundaryConditions();
279 
280  return trho;
281 }
282 
283 
285 {
287  (
288  new volScalarField
289  (
290  IOobject
291  (
292  type() + ':' + mu0_.name_,
293  owner().time().timeName(),
294  owner().regionMesh(),
297  ),
298  owner().regionMesh(),
300  zeroGradientFvPatchScalarField::typeName
301  )
302  );
303 
304  tmu().internalField() = this->mu(0, 0);
305  tmu().correctBoundaryConditions();
306 
307  return tmu;
308 }
309 
310 
312 {
313  tmp<volScalarField> tsigma
314  (
315  new volScalarField
316  (
317  IOobject
318  (
319  type() + ':' + sigma0_.name_,
320  owner().time().timeName(),
321  owner().regionMesh(),
324  ),
325  owner().regionMesh(),
326  dimensionedScalar("0", dimMass/sqr(dimTime), 0.0),
327  zeroGradientFvPatchScalarField::typeName
328  )
329  );
330 
331  tsigma().internalField() = this->sigma(0, 0);
332  tsigma().correctBoundaryConditions();
333 
334  return tsigma;
335 }
336 
337 
339 {
341  (
342  new volScalarField
343  (
344  IOobject
345  (
346  type() + ':' + Cp0_.name_,
347  owner().time().timeName(),
348  owner().regionMesh(),
351  ),
352  owner().regionMesh(),
354  zeroGradientFvPatchScalarField::typeName
355  )
356  );
357 
358  tCp().internalField() = this->Cp(0, 0);
359  tCp().correctBoundaryConditions();
360 
361  return tCp;
362 }
363 
364 
366 {
367  tmp<volScalarField> tkappa
368  (
369  new volScalarField
370  (
371  IOobject
372  (
373  type() + ':' + kappa0_.name_,
374  owner().time().timeName(),
375  owner().regionMesh(),
378  ),
379  owner().regionMesh(),
381  zeroGradientFvPatchScalarField::typeName
382  )
383  );
384 
385  tkappa().internalField() = this->kappa(0, 0);
386  tkappa().correctBoundaryConditions();
387 
388  return tkappa;
389 }
390 
391 
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393 
394 } // End namespace surfaceFilmModels
395 } // End namespace regionModels
396 } // End namespace Foam
397 
398 // ************************************************************************* //
Foam::dimPressure
const dimensionSet dimPressure
constantFilmThermo.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::name_
word name_
Definition: constantFilmThermo.H:62
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData
Definition: constantFilmThermo.H:59
Foam::regionModels::surfaceFilmModels::constantFilmThermo::mu
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
Definition: constantFilmThermo.C:284
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::regionModels::surfaceFilmModels::constantFilmThermo::rho
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
Definition: constantFilmThermo.C:257
Foam::dimDensity
const dimensionSet dimDensity
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Cp
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: constantFilmThermo.C:338
Foam::subModelBase::coeffDict_
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:79
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
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
Definition: filmSubModelBaseI.H:37
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Tb0_
thermoData Tb0_
Boiling temperature [K].
Definition: constantFilmThermo.H:120
Foam::regionModels::surfaceFilmModels::constantFilmThermo::W0_
thermoData W0_
Molecular weight [kg/kmol].
Definition: constantFilmThermo.H:117
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::value_
scalar value_
Definition: constantFilmThermo.H:63
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::set_
bool set_
Definition: constantFilmThermo.H:64
Foam::regionModels::surfaceFilmModels::constantFilmThermo::pv
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
Definition: constantFilmThermo.C:218
Foam::regionModels::surfaceFilmModels::constantFilmThermo::D
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
Definition: constantFilmThermo.C:186
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::regionModels::surfaceFilmModels::constantFilmThermo::kappa0_
thermoData kappa0_
Thermal conductivity [W/m/K].
Definition: constantFilmThermo.H:105
Foam::regionModels::surfaceFilmModels::constantFilmThermo::sigma0_
thermoData sigma0_
Surface tension [kg/s2].
Definition: constantFilmThermo.H:99
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::regionModels::surfaceFilmModels::constantFilmThermo::name_
word name_
Specie name.
Definition: constantFilmThermo.H:90
Foam::regionModels::surfaceFilmModels::constantFilmThermo::mu0_
thermoData mu0_
Dynamic viscosity [Pa.s].
Definition: constantFilmThermo.H:96
Foam::regionModels::surfaceFilmModels::constantFilmThermo::W
virtual scalar W() const
Return molecular weight [kg/kmol].
Definition: constantFilmThermo.C:233
Foam::dimPower
const dimensionSet dimPower
dict
dictionary dict
Definition: searchingEngine.H:14
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
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::regionModels::surfaceFilmModels::constantFilmThermo::rho0_
thermoData rho0_
Density [kg/m3].
Definition: constantFilmThermo.H:93
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::regionModels::surfaceFilmModels::filmThermoModel
Definition: filmThermoModel.H:54
Foam::regionModels::surfaceFilmModels::constantFilmThermo::constantFilmThermo
constantFilmThermo(const constantFilmThermo &)
Disallow default bitwise copy construct.
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::regionModels::surfaceFilmModels::constantFilmThermo::sigma
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
Definition: constantFilmThermo.C:311
filmThermoModel
Base class for film thermo models.
Foam::regionModels::surfaceFilmModels::constantFilmThermo::init
void init(thermoData &td)
Initialise thermoData object.
Definition: constantFilmThermo.C:50
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Cp0_
thermoData Cp0_
Specific heat capacity [J/kg/K].
Definition: constantFilmThermo.H:102
constantFilmThermo
Constant thermo model.
Foam::regionModels::surfaceFilmModels::constantFilmThermo::hl
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
Definition: constantFilmThermo.C:202
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Tb
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
Definition: constantFilmThermo.C:245
timeName
word timeName
Definition: getTimeIndex.H:3
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::regionModels::surfaceFilmModels::constantFilmThermo::kappa
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: constantFilmThermo.C:365
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::regionModels::surfaceFilmModels::constantFilmThermo::~constantFilmThermo
virtual ~constantFilmThermo()
Destructor.
Definition: constantFilmThermo.C:93
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::surfaceFilmModels::constantFilmThermo::name
virtual const word & name() const
Return the specie name.
Definition: constantFilmThermo.C:99