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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "constantFilmThermo.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace regionModels
37{
38namespace surfaceFilmModels
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
44
46(
50);
51
52
53void constantFilmThermo::init(thermoData& td)
54{
55 if (coeffDict_.readIfPresent(td.name_, td.value_))
56 {
57 td.set_ = true;
58 }
59}
60
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
65(
67 const dictionary& dict
68)
69:
70 filmThermoModel(typeName, film, dict),
71 name_(coeffDict_.lookup("specie")),
72 rho0_("rho0"),
73 mu0_("mu0"),
74 sigma0_("sigma0"),
75 Cp0_("Cp0"),
76 kappa0_("kappa0"),
77 D0_("D0"),
78 hl0_("hl0"),
79 pv0_("pv0"),
80 W0_("W0"),
81 Tb0_("Tb0")
82{
83 init(rho0_);
84 init(mu0_);
85 init(sigma0_);
86 init(Cp0_);
87 init(kappa0_);
88 init(D0_);
89 init(hl0_);
90 init(pv0_);
91 init(W0_);
92 init(Tb0_);
93}
94
95
96// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97
99{}
100
101
102// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
103
105{
106 return name_;
107}
108
109
111(
112 const scalar p,
113 const scalar T
114) const
115{
116 if (!rho0_.set_)
117 {
118 coeffDict_.readEntry(rho0_.name_, rho0_.value_);
119 rho0_.set_ = true;
120 }
121
122 return rho0_.value_;
123}
124
125
127(
128 const scalar p,
129 const scalar T
130) const
131{
132 if (!mu0_.set_)
133 {
134 coeffDict_.readEntry(mu0_.name_, mu0_.value_);
135 mu0_.set_ = true;
136 }
137
138 return mu0_.value_;
139}
140
141
143(
144 const scalar p,
145 const scalar T
146) const
147{
148 if (!sigma0_.set_)
149 {
150 coeffDict_.readEntry(sigma0_.name_, sigma0_.value_);
151 sigma0_.set_ = true;
152 }
153
154 return sigma0_.value_;
155}
156
157
159(
160 const scalar p,
161 const scalar T
162) const
163{
164 if (!Cp0_.set_)
165 {
166 coeffDict_.readEntry(Cp0_.name_, Cp0_.value_);
167 Cp0_.set_ = true;
168 }
169
170 return Cp0_.value_;
171}
172
173
175(
176 const scalar p,
177 const scalar T
178) const
179{
180 if (!kappa0_.set_)
181 {
182 coeffDict_.readEntry(kappa0_.name_, kappa0_.value_);
183 kappa0_.set_ = true;
184 }
185
186 return kappa0_.value_;
187}
188
189
191(
192 const scalar p,
193 const scalar T
194) const
195{
196 if (!D0_.set_)
197 {
199 D0_.set_ = true;
200 }
201
202 return D0_.value_;
203}
204
205
207(
208 const scalar p,
209 const scalar T
210) const
211{
212 if (!hl0_.set_)
213 {
214 coeffDict_.readEntry(hl0_.name_, hl0_.value_);
215 hl0_.set_ = true;
216 }
217
218 return hl0_.value_;
219}
220
221
223(
224 const scalar p,
225 const scalar T
226) const
227{
228 if (!pv0_.set_)
229 {
230 coeffDict_.readEntry(pv0_.name_, pv0_.value_);
231 pv0_.set_ = true;
232 }
233
234 return pv0_.value_;
235}
236
237
239{
240 if (!W0_.set_)
241 {
243 W0_.set_ = true;
244 }
245
246 return W0_.value_;
247}
248
249
250scalar constantFilmThermo::Tb(const scalar p) const
251{
252 if (!Tb0_.set_)
253 {
254 coeffDict_.readEntry(Tb0_.name_, Tb0_.value_);
255 Tb0_.set_ = true;
256 }
257
258 return Tb0_.value_;
259}
260
261
263{
265 (
267 (
269 (
270 type() + ':' + rho0_.name_,
271 film().time().timeName(),
272 film().regionMesh(),
275 ),
276 film().regionMesh(),
278 extrapolatedCalculatedFvPatchScalarField::typeName
279 )
280 );
281
282 trho.ref().primitiveFieldRef() = this->rho(0, 0);
283 trho.ref().correctBoundaryConditions();
284
285 return trho;
286}
287
288
290{
292 (
294 (
296 (
297 type() + ':' + mu0_.name_,
298 film().time().timeName(),
299 film().regionMesh(),
302 ),
303 film().regionMesh(),
305 extrapolatedCalculatedFvPatchScalarField::typeName
306 )
307 );
308
309 tmu.ref().primitiveFieldRef() = this->mu(0, 0);
310 tmu.ref().correctBoundaryConditions();
311
312 return tmu;
313}
314
315
317{
319 (
321 (
323 (
324 type() + ':' + sigma0_.name_,
325 film().time().timeName(),
326 film().regionMesh(),
329 ),
330 film().regionMesh(),
332 extrapolatedCalculatedFvPatchScalarField::typeName
333 )
334 );
335
336 tsigma.ref().primitiveFieldRef() = this->sigma(0, 0);
337 tsigma.ref().correctBoundaryConditions();
338
339 return tsigma;
340}
341
342
344{
346 (
348 (
350 (
351 type() + ':' + Cp0_.name_,
352 film().time().timeName(),
353 film().regionMesh(),
356 ),
357 film().regionMesh(),
359 extrapolatedCalculatedFvPatchScalarField::typeName
360 )
361 );
362
363 tCp.ref().primitiveFieldRef() = this->Cp(0, 0);
364 tCp.ref().correctBoundaryConditions();
365
366 return tCp;
367}
368
369
371{
373 (
375 (
377 (
378 type() + ':' + kappa0_.name_,
379 film().time().timeName(),
380 film().regionMesh(),
383 ),
384 film().regionMesh(),
386 extrapolatedCalculatedFvPatchScalarField::typeName
387 )
388 );
389
390 tkappa.ref().primitiveFieldRef() = this->kappa(0, 0);
391 tkappa.ref().correctBoundaryConditions();
392
393 return tkappa;
394}
395
396
397// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398
399} // End namespace surfaceFilmModels
400} // End namespace regionModels
401} // End namespace Foam
402
403// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
Lookup type of boundary radiation properties.
Definition: lookup.H:66
virtual scalar W() const
Return molecular weight [kg/kmol].
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
virtual const word & name() const
Return the specie name.
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:82
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimPower
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict
tmp< volScalarField > trho
const dimensionedScalar & D