Loading...
Searching...
No Matches
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 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "constantFilmThermo.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace regionModels
38{
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47(
51);
52
53
54void constantFilmThermo::init(thermoData& td)
55{
56 if (coeffDict_.readIfPresent(td.name_, td.value_))
57 {
58 td.set_ = true;
59 }
60}
61
62
63// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64
65constantFilmThermo::constantFilmThermo
66(
67 surfaceFilmRegionModel& film,
68 const dictionary& dict
69)
70:
71 filmThermoModel(typeName, film, dict),
72 name_(coeffDict_.lookup("specie")),
73 rho0_("rho0"),
74 mu0_("mu0"),
75 sigma0_("sigma0"),
76 Cp0_("Cp0"),
77 kappa0_("kappa0"),
78 D0_("D0"),
79 hl0_("hl0"),
80 pv0_("pv0"),
81 W0_("W0"),
82 Tb0_("Tb0")
83{
84 init(rho0_);
85 init(mu0_);
86 init(sigma0_);
87 init(Cp0_);
88 init(kappa0_);
89 init(D0_);
90 init(hl0_);
91 init(pv0_);
92 init(W0_);
93 init(Tb0_);
94}
95
96
97// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102
103// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105const word& constantFilmThermo::name() const
106{
107 return name_;
108}
109
110
112(
113 const scalar p,
114 const scalar T
115) const
116{
117 if (!rho0_.set_)
118 {
119 coeffDict_.readEntry(rho0_.name_, rho0_.value_);
120 rho0_.set_ = true;
121 }
122
123 return rho0_.value_;
124}
125
126
128(
129 const scalar p,
130 const scalar T
131) const
132{
133 if (!mu0_.set_)
134 {
135 coeffDict_.readEntry(mu0_.name_, mu0_.value_);
136 mu0_.set_ = true;
137 }
138
139 return mu0_.value_;
140}
141
142
144(
145 const scalar p,
146 const scalar T
147) const
148{
149 if (!sigma0_.set_)
150 {
151 coeffDict_.readEntry(sigma0_.name_, sigma0_.value_);
152 sigma0_.set_ = true;
153 }
154
155 return sigma0_.value_;
156}
157
158
160(
161 const scalar p,
162 const scalar T
163) const
164{
165 if (!Cp0_.set_)
166 {
167 coeffDict_.readEntry(Cp0_.name_, Cp0_.value_);
168 Cp0_.set_ = true;
169 }
170
171 return Cp0_.value_;
172}
173
174
176(
177 const scalar p,
178 const scalar T
179) const
180{
181 if (!kappa0_.set_)
182 {
183 coeffDict_.readEntry(kappa0_.name_, kappa0_.value_);
184 kappa0_.set_ = true;
185 }
186
187 return kappa0_.value_;
188}
189
190
192(
193 const scalar p,
194 const scalar T
195) const
196{
197 if (!D0_.set_)
198 {
200 D0_.set_ = true;
201 }
202
203 return D0_.value_;
204}
205
206
208(
209 const scalar p,
210 const scalar T
211) const
212{
213 if (!hl0_.set_)
214 {
215 coeffDict_.readEntry(hl0_.name_, hl0_.value_);
216 hl0_.set_ = true;
217 }
218
219 return hl0_.value_;
220}
221
222
224(
225 const scalar p,
226 const scalar T
227) const
228{
229 if (!pv0_.set_)
230 {
231 coeffDict_.readEntry(pv0_.name_, pv0_.value_);
232 pv0_.set_ = true;
233 }
234
235 return pv0_.value_;
236}
237
238
239scalar constantFilmThermo::W() const
240{
241 if (!W0_.set_)
242 {
244 W0_.set_ = true;
245 }
246
247 return W0_.value_;
248}
249
250
251scalar constantFilmThermo::Tb(const scalar p) const
252{
253 if (!Tb0_.set_)
254 {
255 coeffDict_.readEntry(Tb0_.name_, Tb0_.value_);
256 Tb0_.set_ = true;
257 }
258
259 return Tb0_.value_;
260}
261
262
264{
266 (
269 film().regionMesh(),
272 );
273
274 trho.ref().primitiveFieldRef() = this->rho(0, 0);
275 trho.ref().correctBoundaryConditions();
276
277 return trho;
278}
279
280
282{
283 auto tmu = volScalarField::New
284 (
287 film().regionMesh(),
289 extrapolatedCalculatedFvPatchScalarField::typeName
290 );
291
292 tmu.ref().primitiveFieldRef() = this->mu(0, 0);
293 tmu.ref().correctBoundaryConditions();
294
295 return tmu;
296}
297
298
300{
301 auto tsigma = volScalarField::New
302 (
303 IOobject::scopedName(type(), sigma0_.name_),
305 film().regionMesh(),
307 extrapolatedCalculatedFvPatchScalarField::typeName
308 );
309
310 tsigma.ref().primitiveFieldRef() = this->sigma(0, 0);
311 tsigma.ref().correctBoundaryConditions();
312
313 return tsigma;
314}
315
316
318{
320 (
323 film().regionMesh(),
325 extrapolatedCalculatedFvPatchScalarField::typeName
326 );
327
328 tCp.ref().primitiveFieldRef() = this->Cp(0, 0);
329 tCp.ref().correctBoundaryConditions();
330
331 return tCp;
332}
333
334
336{
337 auto tkappa = volScalarField::New
338 (
339 IOobject::scopedName(type(), kappa0_.name_),
341 film().regionMesh(),
343 extrapolatedCalculatedFvPatchScalarField::typeName
344 );
345
346 tkappa.ref().primitiveFieldRef() = this->kappa(0, 0);
347 tkappa.ref().correctBoundaryConditions();
348
349 return tkappa;
350}
351
352
353// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354
355} // End namespace surfaceFilmModels
356} // End namespace regionModels
357} // End namespace Foam
358
359// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Lookup type of boundary radiation properties.
Definition lookup.H:60
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 scalar Tb(const scalar p) const
Return boiling temperature [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].
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const tmp< volScalarField > & tCp
Definition EEqn.H:4
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimPower
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
tmp< volScalarField > trho