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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 namespace areaSurfaceFilmModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  scalarField X(thermo_.size(), 1);
50 
51  forAll(rho_, faceI)
52  {
53  rho_[faceI] = thermo_.rho(pRef_, Tf_[faceI], X);
54  mu_[faceI] = thermo_.mu(pRef_, Tf_[faceI], X);
55  sigma_[faceI] = thermo_.sigma(pRef_, Tf_[faceI], X);
56  Cp_[faceI] = thermo_.Cp(pRef_, Tf_[faceI], X);
57  }
58 
59  forAll(regionMesh().boundary(), patchI)
60  {
61  const scalarField& patchTf = Tf_.boundaryFieldRef()[patchI];
62 
63  scalarField& patchRho = rho_.boundaryFieldRef()[patchI];
64  scalarField& patchmu = mu_.boundaryFieldRef()[patchI];
65  scalarField& patchsigma = sigma_.boundaryFieldRef()[patchI];
66  scalarField& patchCp = Cp_.boundaryFieldRef()[patchI];
67 
68  forAll(patchRho, edgeI)
69  {
70  patchRho[edgeI] = thermo_.rho(pRef_, patchTf[edgeI], X);
71  patchmu[edgeI] = thermo_.mu(pRef_, patchTf[edgeI], X);
72  patchsigma[edgeI] = thermo_.sigma(pRef_, patchTf[edgeI], X);
73  patchCp[edgeI] = thermo_.Cp(pRef_, patchTf[edgeI], X);
74  }
75  }
76 
77  //Initialize pf_
79 }
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const word& modelType,
86  const fvMesh& mesh,
87  const dictionary& dict
88 )
89 :
90  liquidFilmBase(modelType, mesh, dict),
91  thermo_(dict.subDict("thermo")),
92  rho_
93  (
94  IOobject
95  (
96  "rhof",
97  regionMesh().time().timeName(),
98  regionMesh().thisDb(),
99  IOobject::NO_READ,
100  IOobject::AUTO_WRITE
101  ),
102  regionMesh(),
104  ),
105  mu_
106  (
107  IOobject
108  (
109  "muf",
110  regionMesh().time().timeName(),
111  regionMesh().thisDb(),
112  IOobject::NO_READ,
113  IOobject::NO_WRITE
114  ),
115  regionMesh(),
117  ),
118  Tf_
119  (
120  IOobject
121  (
122  "Tf_" + regionName_,
123  regionMesh().time().timeName(),
124  regionMesh().thisDb(),
125  IOobject::READ_IF_PRESENT,
126  IOobject::AUTO_WRITE
127  ),
128  regionMesh(),
130  ),
131  Cp_
132  (
133  IOobject
134  (
135  "Cp_" + regionName_,
136  regionMesh().time().timeName(),
137  regionMesh().thisDb(),
138  IOobject::NO_READ,
139  IOobject::NO_WRITE
140  ),
141  regionMesh(),
143  ),
144  sigma_
145  (
146  IOobject
147  (
148  "sigmaf",
149  regionMesh().time().timeName(),
150  regionMesh().thisDb(),
151  IOobject::NO_READ,
152  IOobject::NO_WRITE
153  ),
154  regionMesh(),
156  ),
157  hRho_
158  (
159  IOobject
160  (
161  h_.name() + "*" + rho_.name(),
162  regionMesh().time().timeName(),
163  regionMesh().thisDb(),
164  IOobject::NO_READ,
165  IOobject::NO_WRITE
166  ),
167  regionMesh(),
168  dimensionedScalar(h_.dimensions()*rho_.dimensions(), Zero)
169  ),
170  rhoSp_
171  (
172  IOobject
173  (
174  "rhoSp",
175  regionMesh().time().timeName(),
176  regionMesh().thisDb(),
177  IOobject::NO_READ,
178  IOobject::NO_WRITE
179  ),
180  regionMesh(),
182  ),
183  USp_
184  (
185  IOobject
186  (
187  "USp",
188  regionMesh().time().timeName(),
189  regionMesh().thisDb(),
190  IOobject::NO_READ,
191  IOobject::NO_WRITE
192  ),
193  regionMesh(),
195  ),
196  pnSp_
197  (
198  IOobject
199  (
200  "pnSp",
201  regionMesh().time().timeName(),
202  regionMesh().thisDb(),
203  IOobject::NO_READ,
204  IOobject::NO_WRITE
205  ),
206  regionMesh(),
208  ),
209  cloudMassTrans_
210  (
211  IOobject
212  (
213  "cloudMassTrans",
214  primaryMesh().time().timeName(),
215  primaryMesh().thisDb(),
216  IOobject::NO_READ,
217  IOobject::NO_WRITE
218  ),
219  primaryMesh(),
221  ),
222  cloudDiameterTrans_
223  (
224  IOobject
225  (
226  "cloudDiameterTrans",
227  primaryMesh().time().timeName(),
228  primaryMesh().thisDb(),
229  IOobject::NO_READ,
230  IOobject::NO_WRITE
231  ),
232  primaryMesh(),
234  ),
235  turbulence_(filmTurbulenceModel::New(*this, dict)),
236  availableMass_(regionMesh().faces().size(), Zero),
237  injection_(*this, dict),
238  forces_(*this, dict)
239 {
240  if (dict.readIfPresent("T0", Tref_))
241  {
242  Info<< "Initialised film temperature to T0" << endl;
243 
245  }
247 }
248 
249 
250 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
253 {
254  return cloudMassTrans_;
255 }
256 
259 {
260  return cloudDiameterTrans_;
261 }
262 
263 
265 {
267 
268 
271 
272  const scalar deltaT = primaryMesh().time().deltaTValue();
273  const scalarField rAreaDeltaT(scalar(1)/deltaT/regionMesh().S().field());
274 
275  // Map the total mass, mom [kg.m/s] and pnSource from particles
276 
278 
280 
282 
283 
284  // Calculate rate per area
285  rhoSp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
286  USp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
287  pnSp_.primitiveFieldRef() *= rAreaDeltaT/rho_;
289  rhoSp_.relax();
290  pnSp_.relax();
291  USp_.relax();
292 }
293 
294 
296 {
297  availableMass_ = (h() - h0_)*rho()*regionMesh().S();
300 }
301 
302 
304 {
305  Info<< "\nSurface film: " << type() << " on patch: ";
306 
307  for (const label patchi : this->primaryPatchIDs())
308  {
309  Info<< ' ' << patchi;
310  }
311  Info<< endl;
312 
314 
315  Info<< indent << "min/max(mag(Uf)) = "
316  << gMinMaxMag(Uf_.field()) << nl
317  << indent << "min/max(delta) = "
318  << gMinMax(h_.field()) << nl
319  << indent << "coverage = "
320  << gSum(alpha()().field()*mag(sf.field()))/gSumMag(sf.field()) << nl
321  << indent << "total mass = "
322  << gSum(availableMass_) << nl;
323 
325 }
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 } // End namespace areaSurfaceFilmModels
330 } // End namespace regionModels
331 } // End namespace Foam
332 
333 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
faceListList boundary
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
dictionary dict
rDeltaTY field()
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
areaScalarField mu_
Dynamic viscosity [Pa.s].
const areaScalarField & rho() const noexcept
Access const reference rho.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar h0_
Smallest film thickness.
label size() const
Return the number of liquids in the mixture.
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:1133
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
virtual void info(Ostream &os)
Provide some info.
defineTypeNameAndDebug(kinematicThinFilm, 0)
Macros for easy insertion into run-time selection tables.
volScalarField pnSource_
Normal pressure by particles.
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
areaScalarField sigma_
Surface tension [m/s^2].
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
const areaScalarField & h() const noexcept
Access const reference h.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
const faMesh & regionMesh() const
Return the region mesh database.
volScalarField cloudMassTrans_
Film mass for transfer to cloud.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
const dimensionSet dimPressure
areaScalarField Cp_
Film heat capacity [J/K].
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...
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
liquidFilmModel(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
scalarMinMax gMinMaxMag(const FieldField< Field, Type > &f)
const dimensionSet dimEnergy
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
const dimensionSet dimDensity
const labelList & primaryPatchIDs() const
List of patch IDs on the primary region coupled to this region.
const Field< Type > & field() const noexcept
Return const-reference to the field values.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
areaScalarField pnSp_
Normal pressure by particles.
void relax(const scalar alpha)
Relax field (for steady-state solution).
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
scalarField availableMass_
Available mass for transfer via sub-models.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
tmp< areaScalarField > alpha() const
Wet indicator using h0.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity