constantFilmThermo.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 "constantFilmThermo.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(constantFilmThermo, 0);
42 
44 (
45  filmThermoModel,
46  constantFilmThermo,
47  dictionary
48 );
49 
50 
51 void constantFilmThermo::init(thermoData& td)
52 {
53  if (coeffDict_.readIfPresent(td.name_, td.value_))
54  {
55  td.set_ = true;
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
62 constantFilmThermo::constantFilmThermo
63 (
65  const dictionary& dict
66 )
67 :
68  filmThermoModel(typeName, film, dict),
69  name_(coeffDict_.lookup("specie")),
70  rho0_("rho0"),
71  mu0_("mu0"),
72  sigma0_("sigma0"),
73  Cp0_("Cp0"),
74  kappa0_("kappa0"),
75  D0_("D0"),
76  hl0_("hl0"),
77  pv0_("pv0"),
78  W0_("W0"),
79  Tb0_("Tb0")
80 {
81  init(rho0_);
82  init(mu0_);
83  init(sigma0_);
84  init(Cp0_);
85  init(kappa0_);
86  init(D0_);
87  init(hl0_);
88  init(pv0_);
89  init(W0_);
90  init(Tb0_);
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
97 {}
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
103 {
104  return name_;
105 }
106 
107 
109 (
110  const scalar p,
111  const scalar T
112 ) const
113 {
114  if (!rho0_.set_)
115  {
116  coeffDict_.lookup(rho0_.name_) >> rho0_.value_;
117  rho0_.set_ = true;
118  }
119 
120  return rho0_.value_;
121 }
122 
123 
125 (
126  const scalar p,
127  const scalar T
128 ) const
129 {
130  if (!mu0_.set_)
131  {
132  coeffDict_.lookup(mu0_.name_) >> mu0_.value_;
133  mu0_.set_ = true;
134  }
135 
136  return mu0_.value_;
137 }
138 
139 
141 (
142  const scalar p,
143  const scalar T
144 ) const
145 {
146  if (!sigma0_.set_)
147  {
148  coeffDict_.lookup(sigma0_.name_) >> sigma0_.value_;
149  sigma0_.set_ = true;
150  }
151 
152  return sigma0_.value_;
153 }
154 
155 
157 (
158  const scalar p,
159  const scalar T
160 ) const
161 {
162  if (!Cp0_.set_)
163  {
164  coeffDict_.lookup(Cp0_.name_) >> Cp0_.value_;
165  Cp0_.set_ = true;
166  }
167 
168  return Cp0_.value_;
169 }
170 
171 
173 (
174  const scalar p,
175  const scalar T
176 ) const
177 {
178  if (!kappa0_.set_)
179  {
180  coeffDict_.lookup(kappa0_.name_) >> kappa0_.value_;
181  kappa0_.set_ = true;
182  }
183 
184  return kappa0_.value_;
185 }
186 
187 
189 (
190  const scalar p,
191  const scalar T
192 ) const
193 {
194  if (!D0_.set_)
195  {
196  coeffDict_.lookup(D0_.name_) >> D0_.value_;
197  D0_.set_ = true;
198  }
199 
200  return D0_.value_;
201 }
202 
203 
205 (
206  const scalar p,
207  const scalar T
208 ) const
209 {
210  if (!hl0_.set_)
211  {
212  coeffDict_.lookup(hl0_.name_) >> hl0_.value_;
213  hl0_.set_ = true;
214  }
215 
216  return hl0_.value_;
217 }
218 
219 
221 (
222  const scalar p,
223  const scalar T
224 ) const
225 {
226  if (!pv0_.set_)
227  {
228  coeffDict_.lookup(pv0_.name_) >> pv0_.value_;
229  pv0_.set_ = true;
230  }
231 
232  return pv0_.value_;
233 }
234 
235 
236 scalar constantFilmThermo::W() const
237 {
238  if (!W0_.set_)
239  {
240  coeffDict_.lookup(W0_.name_) >> W0_.value_;
241  W0_.set_ = true;
242  }
243 
244  return W0_.value_;
245 }
246 
247 
248 scalar constantFilmThermo::Tb(const scalar p) const
249 {
250  if (!Tb0_.set_)
251  {
252  coeffDict_.lookup(Tb0_.name_) >> Tb0_.value_;
253  Tb0_.set_ = true;
254  }
255 
256  return Tb0_.value_;
257 }
258 
259 
261 {
263  (
264  new volScalarField
265  (
266  IOobject
267  (
268  type() + ':' + rho0_.name_,
269  film().time().timeName(),
270  film().regionMesh(),
273  ),
274  film().regionMesh(),
275  dimensionedScalar("0", dimDensity, 0.0),
276  extrapolatedCalculatedFvPatchScalarField::typeName
277  )
278  );
279 
280  trho.ref().primitiveFieldRef() = this->rho(0, 0);
281  trho.ref().correctBoundaryConditions();
282 
283  return trho;
284 }
285 
286 
288 {
290  (
291  new volScalarField
292  (
293  IOobject
294  (
295  type() + ':' + mu0_.name_,
296  film().time().timeName(),
297  film().regionMesh(),
300  ),
301  film().regionMesh(),
303  extrapolatedCalculatedFvPatchScalarField::typeName
304  )
305  );
306 
307  tmu.ref().primitiveFieldRef() = this->mu(0, 0);
308  tmu.ref().correctBoundaryConditions();
309 
310  return tmu;
311 }
312 
313 
315 {
316  tmp<volScalarField> tsigma
317  (
318  new volScalarField
319  (
320  IOobject
321  (
322  type() + ':' + sigma0_.name_,
323  film().time().timeName(),
324  film().regionMesh(),
327  ),
328  film().regionMesh(),
329  dimensionedScalar("0", dimMass/sqr(dimTime), 0.0),
330  extrapolatedCalculatedFvPatchScalarField::typeName
331  )
332  );
333 
334  tsigma.ref().primitiveFieldRef() = this->sigma(0, 0);
335  tsigma.ref().correctBoundaryConditions();
336 
337  return tsigma;
338 }
339 
340 
342 {
344  (
345  new volScalarField
346  (
347  IOobject
348  (
349  type() + ':' + Cp0_.name_,
350  film().time().timeName(),
351  film().regionMesh(),
354  ),
355  film().regionMesh(),
357  extrapolatedCalculatedFvPatchScalarField::typeName
358  )
359  );
360 
361  tCp.ref().primitiveFieldRef() = this->Cp(0, 0);
362  tCp.ref().correctBoundaryConditions();
363 
364  return tCp;
365 }
366 
367 
369 {
370  tmp<volScalarField> tkappa
371  (
372  new volScalarField
373  (
374  IOobject
375  (
376  type() + ':' + kappa0_.name_,
377  film().time().timeName(),
378  film().regionMesh(),
381  ),
382  film().regionMesh(),
384  extrapolatedCalculatedFvPatchScalarField::typeName
385  )
386  );
387 
388  tkappa.ref().primitiveFieldRef() = this->kappa(0, 0);
389  tkappa.ref().correctBoundaryConditions();
390 
391  return tkappa;
392 }
393 
394 
395 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396 
397 } // End namespace surfaceFilmModels
398 } // End namespace regionModels
399 } // End namespace Foam
400 
401 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volScalarField > trho
virtual const word & name() const
Return the specie name.
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:79
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
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
virtual scalar W() const
Return molecular weight [kg/kmol].
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimPressure
const dimensionSet dimPower
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
const dimensionSet dimEnergy
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
const dimensionSet dimDensity
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576