homogeneousCondensation.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) 2024 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 
27 #include "fundamentalConstants.H"
28 #include "multicomponentThermo.H"
30 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
41 }
42 }
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::fv::homogeneousCondensation::readCoeffs()
47 {
48  saturationModel_.reset
49  (
50  saturationPressureModel::New("pSat", coeffs()).ptr()
51  );
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
58 (
59  const word& name,
60  const word& modelType,
61  const fvMesh& mesh,
62  const dictionary& dict
63 )
64 :
66  (
67  name,
68  modelType,
69  mesh,
70  dict,
71  {false, false},
72  {true, false}
73  ),
74  fluid_
75  (
76  mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)
77  ),
78  interface_
79  (
80  fluid_.phases()[phaseNames().first()],
81  fluid_.phases()[phaseNames().second()]
82  ),
83  d_
84  (
85  IOobject
86  (
87  typedName(IOobject::groupName("d", interface_.name())),
88  mesh.time().name(),
89  mesh,
92  ),
93  mesh,
95  ),
96  mDotByAlphaGas_
97  (
98  IOobject
99  (
100  typedName(IOobject::groupName("mDotByAlpha", interface_.name())),
101  mesh.time().name(),
102  mesh,
105  ),
106  mesh,
108  ),
109  saturationModel_(nullptr)
110 {
111  readCoeffs();
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
119 {
120  return d_;
121 }
122 
123 
126 {
127  const volScalarField::Internal& alphaGas =
128  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
129 
130  const multicomponentThermo& thermoGas = specieThermos().first();
131 
132  const volScalarField& p = this->p();
133  const volScalarField& T = thermoGas.T();
134 
135  const volScalarField::Internal rhoLiquid
136  (
137  specieThermos().valid().second()
138  ? vfToVif(specieThermos().second().rhoi(specieis().second(), p, T))
139  : vfToVif(thermos().second().rho())
140  );
141 
143 
144  return alphaGas*mDotByAlphaGas_/(rhoLiquid*v);
145 }
146 
147 
150 {
151  const volScalarField::Internal& alphaGas =
152  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
153 
154  return alphaGas*mDotByAlphaGas_;
155 }
156 
157 
159 {
160  #define DebugField(field) \
161  DebugInfo \
162  << name() << ": " \
163  << #field << ' ' << field.dimensions() << " min/avg/max = " \
164  << gMin(field) << '/' << gAverage(field) << '/' << gMax(field) \
165  << nl;
166 
170 
171  const multicomponentThermo& thermoGas = specieThermos().first();
172 
173  const volScalarField& p = this->p();
174  const volScalarField& T = thermoGas.T();
175  DebugField(p);
176  DebugField(T);
177 
178  // Phase molecular masses and densities
179  const volScalarField::Internal WGas(vfToVif(thermoGas.W()));
180  const volScalarField::Internal rhoGas(vfToVif(thermoGas.rho()));
181  const volScalarField::Internal WLiquid
182  (
183  specieThermos().valid().second()
185  (
186  "W",
187  mesh(),
188  specieThermos().second().Wi(specieis().second())
189  )
190  : vfToVif(thermos().second().W())
191  );
192  const volScalarField::Internal rhoLiquid
193  (
194  specieThermos().valid().second()
195  ? vfToVif(specieThermos().second().rhoi(specieis().second(), p, T))
196  : vfToVif(thermos().second().rho())
197  );
198  DebugField(WGas);
199  DebugField(rhoGas);
200  DebugField(WLiquid);
201  DebugField(rhoLiquid);
202 
203  // Surface tension
204  const volScalarField::Internal sigma(interface_.sigma());
205  DebugField(sigma);
206 
207  // Mole fraction of nucleating specie
208  const volScalarField::Internal Xi
209  (
210  thermoGas.Y()[specieis().first()]*WGas/thermoGas.Wi(specieis().first())
211  );
212 
213  // Saturation pressure and concentration
214  const volScalarField::Internal pSat(saturationModel_->pSat(T()));
215  const volScalarField::Internal cSat(pSat/p()*rhoGas/WGas);
216  DebugField(pSat);
217  DebugField(cSat);
218 
219  // Supersaturation of the nucleating specie
220  const volScalarField::Internal S(Xi*p()/pSat);
221  DebugField(S);
222 
223  // Mass, volume and diameter of one molecule in the condensed phase
224  const volScalarField::Internal mMolc(WLiquid/NA);
225  const volScalarField::Internal vMolc(mMolc/rhoLiquid);
226  const volScalarField::Internal dMolc(cbrt(6/pi*vMolc));
227  DebugField(mMolc);
228  DebugField(vMolc);
229  DebugField(dMolc);
230 
231  // Diameter of nuclei
232  d_ = 4*sigma*vMolc/(k*T()*log(max(S, 1 + small)));
233  DebugField(d_);
234 
235  // ?
236  const volScalarField::Internal deltaPhiStar(pi/3*sigma*sqr(d_));
237  DebugField(deltaPhiStar);
238 
239  // Ratio of nucleus volume to molecular volume
240  const volScalarField::Internal iStar(pi/6*pow3(d_)/vMolc);
241  DebugField(iStar);
242 
243  // ?
244  const volScalarField::Internal betaIStar1
245  (
246  sqrt(6*k*T()/mMolc)*sqrt((iStar + 1)/iStar)*sqr(d_/2 + dMolc/2)
247  );
248  DebugField(betaIStar1);
249 
250  // Number-based nucleation rate; i.e., number of nuclei created per second
251  // per unit volume
253  (
254  betaIStar1*sqr(cSat*NA)*exp(-deltaPhiStar/(k*T()))
255  *sqrt(sigma/(k*T()))
256  *2*mMolc/(pi*sqr(d_)*rhoLiquid)
257  );
258  DebugField(J);
259 
260  // Mass transfer rate
261  mDotByAlphaGas_ = J*iStar*mMolc;
262  DebugField(mDotByAlphaGas_);
263 
264  #undef DebugField
265 }
266 
267 
269 (
270  const volScalarField& alpha,
271  const volScalarField& rho,
272  fvMatrix<scalar>& eqn
273 ) const
274 {
275  const label i = index(alphaNames(), eqn.psi().name());
276 
277  if (i != -1)
278  {
279  if (i == 0)
280  {
281  eqn -= fvm::Sp(mDotByAlphaGas_, eqn.psi());
282  }
283  else
284  {
285  eqn += mDot() - correction(fvm::Sp(mDotByAlphaGas_, eqn.psi()));
286  }
287  }
288  else
289  {
291  }
292 }
293 
294 
296 (
297  const volScalarField& alpha,
298  const volScalarField& rho,
299  const volScalarField& Yi,
300  fvMatrix<scalar>& eqn
301 ) const
302 {
303  const label i = index(alphaNames(), eqn.psi().name());
304 
305  if (i == 0 && specieis().first() != -1 && Yi.member() == specie())
306  {
307  eqn -= fvm::Sp(alpha*mDotByAlphaGas_, Yi);
308  }
309  else
310  {
312  }
313 }
314 
315 
317 {
319  {
320  readCoeffs();
321  return true;
322  }
323  else
324  {
325  return false;
326  }
327 }
328 
329 
330 // ************************************************************************* //
label k
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...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:149
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
virtual const volScalarField & T() const =0
Temperature [K].
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
Model for the homogeneous nucleation of liquid droplets out of a gaseous mixture.
virtual void correct()
Correct the fvModel.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
homogeneousCondensation(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual tmp< DimensionedField< scalar, volMesh > > d() const
Return the diameter of nuclei.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the compressible continuity equation to add.
virtual tmp< DimensionedField< scalar, volMesh > > nDot() const
Return the number rate at which nuclei are generated.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:218
Base class for phase change models in which only a single component changes phase....
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Base-class for multi-component thermodynamic properties.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
virtual dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:226
static autoPtr< saturationPressureModel > New(const dictionary &dict)
Select with dictionary.
Base class of the thermophysical property types.
Definition: specie.H:66
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Fundamental dimensioned constants.
#define DebugField(field)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar k
Boltzmann constant.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar NA
Avagadro number.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
dimensionedScalar log(const dimensionedScalar &ds)
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList fv(nPoints)
dictionary dict
volScalarField & p