liquidFilmModel.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) 2020-2021 OpenCFD Ltd.
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 "liquidFilmModel.H"
31#include "gravityMeshObject.H"
32#include "volFields.H"
33
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace regionModels
40{
41namespace areaSurfaceFilmModels
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
52 scalarField X(thermo_.size(), 1);
53
54 forAll (rho_, faceI)
55 {
56 rho_[faceI] = thermo_.rho(pRef_, Tf_[faceI], X);
57 mu_[faceI] = thermo_.mu(pRef_, Tf_[faceI], X);
58 sigma_[faceI] = thermo_.sigma(pRef_, Tf_[faceI], X);
59 Cp_[faceI] = thermo_.Cp(pRef_, Tf_[faceI], X);
60 }
61
62 forAll (regionMesh().boundary(), patchI)
63 {
64 const scalarField& patchTf = Tf_.boundaryFieldRef()[patchI];
65
66 scalarField& patchRho = rho_.boundaryFieldRef()[patchI];
67 scalarField& patchmu = mu_.boundaryFieldRef()[patchI];
68 scalarField& patchsigma = sigma_.boundaryFieldRef()[patchI];
69 scalarField& patchCp = Cp_.boundaryFieldRef()[patchI];
70
71 forAll(patchRho, edgeI)
72 {
73 patchRho[edgeI] = thermo_.rho(pRef_, patchTf[edgeI], X);
74 patchmu[edgeI] = thermo_.mu(pRef_, patchTf[edgeI], X);
75 patchsigma[edgeI] = thermo_.sigma(pRef_, patchTf[edgeI], X);
76 patchCp[edgeI] = thermo_.Cp(pRef_, patchTf[edgeI], X);
77 }
78 }
79
80 //Initialize pf_
82}
83
84// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85
87(
88 const word& modelType,
89 const fvPatch& patch,
90 const dictionary& dict
91)
92:
93 liquidFilmBase(modelType, patch, dict),
94
95 thermo_(dict.subDict("thermo")),
96
97 rho_
98 (
100 (
101 "rhof",
102 primaryMesh().time().timeName(),
103 primaryMesh(),
104 IOobject::NO_READ,
105 IOobject::AUTO_WRITE
106 ),
107 regionMesh(),
109 ),
110 mu_
111 (
113 (
114 "muf",
115 primaryMesh().time().timeName(),
116 primaryMesh(),
117 IOobject::NO_READ,
118 IOobject::NO_WRITE
119 ),
120 regionMesh(),
122 ),
123 Tf_
124 (
126 (
127 "Tf_" + regionName_,
128 primaryMesh().time().timeName(),
129 primaryMesh(),
130 IOobject::READ_IF_PRESENT,
131 IOobject::AUTO_WRITE
132 ),
133 regionMesh(),
135 ),
136 Cp_
137 (
139 (
140 "Cp_" + regionName_,
141 primaryMesh().time().timeName(),
142 primaryMesh(),
143 IOobject::NO_READ,
144 IOobject::NO_WRITE
145 ),
146 regionMesh(),
148 ),
149 sigma_
150 (
152 (
153 "sigmaf",
154 primaryMesh().time().timeName(),
155 primaryMesh(),
156 IOobject::NO_READ,
157 IOobject::NO_WRITE
158 ),
159 regionMesh(),
161 ),
162 hRho_
163 (
165 (
166 h_.name() + "*" + rho_.name(),
167 primaryMesh().time().timeName(),
168 primaryMesh(),
169 IOobject::NO_READ,
170 IOobject::NO_WRITE
171 ),
172 regionMesh(),
173 dimensionedScalar(h_.dimensions()*rho_.dimensions(), Zero)
174 ),
175 rhoSp_
176 (
178 (
179 "rhoSp",
180 primaryMesh().time().timeName(),
181 primaryMesh(),
182 IOobject::NO_READ,
183 IOobject::NO_WRITE
184 ),
185 regionMesh(),
187 ),
188 USp_
189 (
191 (
192 "USp",
193 primaryMesh().time().timeName(),
194 primaryMesh()
195 ),
196 regionMesh(),
198 ),
199 pnSp_
200 (
202 (
203 "pnSp",
204 primaryMesh().time().timeName(),
205 primaryMesh()
206 ),
207 regionMesh(),
209 ),
210 cloudMassTrans_
211 (
213 (
214 "cloudMassTrans",
215 primaryMesh().time().timeName(),
216 primaryMesh(),
217 IOobject::NO_READ,
218 IOobject::NO_WRITE
219 ),
220 primaryMesh(),
222 calculatedFvPatchField<scalar>::typeName
223 ),
224
225 cloudDiameterTrans_
226 (
228 (
229 "cloudDiameterTrans",
230 primaryMesh().time().timeName(),
231 primaryMesh(),
232 IOobject::NO_READ,
233 IOobject::NO_WRITE
234 ),
235 primaryMesh(),
237 calculatedFvPatchField<scalar>::typeName
238 ),
239
240 turbulence_(filmTurbulenceModel::New(*this, dict)),
241
242 availableMass_(regionMesh().faces().size(), Zero),
243
244 injection_(*this, dict),
245
246 forces_(*this, dict)
247{
248
249 if (dict.found("T0"))
250 {
251 Tref_ = dict.get<scalar>("T0");
253 }
255}
256
257
258// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
259
261{
262 return mu_;
263}
264
265
267{
268 return rho_;
269}
270
271
273{
274 return sigma_;
275}
276
277
279{
280 return Tf_;
281}
282
283
285{
286 return Cp_;
287}
288
289
291{
292 return thermo_;
293}
294
295
297{
298 return Tref_;
299}
300
301
303{
304 return cloudMassTrans_;
305}
306
307
309{
310 return cloudDiameterTrans_;
311}
312
313
315{
317
318
319
322
323 const scalar deltaT = primaryMesh().time().deltaTValue();
324 const scalarField rAreaDeltaT(scalar(1)/deltaT/regionMesh().S().field());
325
326 // Map the total mass, mom and pnSource from particles
329 // [kg.m/s]
332
335
336
337 // Calculate rate per area
338 rhoSp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
339 USp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
340 pnSp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
341
342 rhoSp_.relax();
343 pnSp_.relax();
344 USp_.relax();
345}
346
347
349{
350 availableMass_ = (h() - h0_)*rho()*regionMesh().S();
353}
354
355
357{
358 Info<< "\nSurface film: " << type() << " on patch: " << patchID() << endl;
359
361
362 Info<< indent << "min/max(mag(Uf)) = "
363 << gMin(mag(Uf_.field())) << ", "
364 << gMax(mag(Uf_.field())) << nl
365 << indent << "min/max(delta) = "
366 << gMin(h_.field()) << ", " << gMax(h_.field()) << nl
367 << indent << "coverage = "
368 << gSum(alpha()().field()*mag(sf.field()))/gSum(mag(sf.field())) << nl
369 << indent << "total mass = "
370 << gSum(availableMass_) << nl;
371
373}
374
375// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
376
377} // End namespace areaSurfaceFilmModels
378} // End namespace regionModels
379} // End namespace Foam
380
381// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Field< Type > & field() const
Return field.
void relax(const scalar alpha)
Relax field (for steady-state solution).
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:780
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
label size() const
Return the number of liquids in the mixture.
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
tmp< areaScalarField > alpha() const
Wet indicator using h0.
volScalarField pnSource_
Normal pressure by particles.
dimensionedScalar h0_
Smallest numerical thickness.
const areaScalarField & h() const
Access const reference h.
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
areaScalarField pnSp_
Normal pressure by particles.
scalar Tref() const
Access to reference temperature.
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
const areaScalarField & Tf() const
Access const reference Tf.
const areaScalarField & Cp() const
Access const reference Cp.
areaScalarField mu_
Dynamic viscosity [Pa.s].
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
const areaScalarField & mu() const
Access const reference mu.
const areaScalarField & rho() const
Access const reference rho.
const areaScalarField & sigma() const
Access const reference sigma.
volScalarField cloudMassTrans_
Film mass for transfer to cloud.
const liquidMixtureProperties & thermo() const
Access to thermo.
scalarField availableMass_
Available mass for transfer via sub-models.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
label patchID() const
Return patch ID.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const faMesh & regionMesh() const
Return the region mesh database.
tmp< Field< Type > > mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &df) const
Map volume boundary field to surface.
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
rDeltaTY field()
faceListList boundary
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimViscosity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gSum(const FieldField< Field, Type > &f)
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
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333