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-2025 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(const dictionary& dict)
47 {
49 
50  saturationModel_.reset
51  (
52  saturationPressureModel::New("pSat", dict).ptr()
53  );
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const word& name,
62  const word& modelType,
63  const fvMesh& mesh,
64  const dictionary& dict
65 )
66 :
67  phaseChange(name, modelType, mesh, dict, readSpecie(dict, true)),
68  fluid_
69  (
70  mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)
71  ),
72  d_
73  (
74  IOobject
75  (
76  name + ":d",
77  mesh.time().name(),
78  mesh,
79  IOobject::READ_IF_PRESENT,
80  IOobject::AUTO_WRITE
81  ),
82  mesh,
84  ),
85  mDotByAlphaGas_
86  (
87  IOobject
88  (
89  name + ":mDotByAlpha",
90  mesh.time().name(),
91  mesh,
92  IOobject::READ_IF_PRESENT,
93  IOobject::AUTO_WRITE
94  ),
95  mesh,
97  ),
98  saturationModel_(nullptr)
99 {
100  readCoeffs(coeffs(dict));
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
108 {
109  return d_;
110 }
111 
112 
115 {
116  const volScalarField::Internal& alphaGas =
117  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
118 
119  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
120  this->multicomponentThermos(true, false);
121 
122  const multicomponentThermo& thermoGas = multicomponentThermos.first();
123 
124  const volScalarField& p = this->p();
125  const volScalarField& T = thermoGas.T();
126 
127  const volScalarField::Internal rhoLiquid
128  (
129  multicomponentThermos.valid().second()
130  ? vfToVif(multicomponentThermos.second().rhoi(specieis().second(), p, T))
131  : vfToVif(thermos().second().rho())
132  );
133 
135 
136  return alphaGas*mDotByAlphaGas_/(rhoLiquid*v);
137 }
138 
139 
142 {
143  const volScalarField::Internal& alphaGas =
144  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
145 
146  return alphaGas*mDotByAlphaGas_;
147 }
148 
149 
152 {
153  static const dimensionedScalar mDotRootVSmall
154  (
156  rootVSmall
157  );
158 
159  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
160  this->multicomponentThermos(true, false);
161 
162  const multicomponentThermo& thermoGas = multicomponentThermos.first();
163 
164  // Phase molecular masses and densities
165  const volScalarField::Internal WGas(vfToVif(thermoGas.W()));
166  const volScalarField::Internal rhoGas(vfToVif(thermoGas.rho()));
167 
168  // Mole fraction of nucleating specie
169  const volScalarField::Internal Xi
170  (
171  thermoGas.Y()[specieis().first()]*WGas/thermoGas.Wi(specieis().first())
172  );
173 
174  return Xi*rhoGas/max(mDotByAlphaGas_, mDotRootVSmall);
175 }
176 
177 
179 {
180  #define DebugField(field) \
181  DebugInfo \
182  << name() << ": " \
183  << #field << ' ' << field.dimensions() << " min/avg/max = " \
184  << gMin(field) << '/' << gAverage(field) << '/' << gMax(field) \
185  << nl;
186 
190 
191  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
192  this->multicomponentThermos(true, false);
193 
194  const multicomponentThermo& thermoGas = multicomponentThermos.first();
195 
196  const volScalarField& p = this->p();
197  const volScalarField& T = thermoGas.T();
198  DebugField(p);
199  DebugField(T);
200 
201  // Phase molecular masses and densities
202  const volScalarField::Internal WGas(vfToVif(thermoGas.W()));
203  const volScalarField::Internal rhoGas(vfToVif(thermoGas.rho()));
204  const volScalarField::Internal WLiquid
205  (
206  multicomponentThermos.valid().second()
208  (
209  "W",
210  mesh(),
211  multicomponentThermos.second().Wi(specieis().second())
212  )
213  : vfToVif(thermos().second().W())
214  );
215  const volScalarField::Internal rhoLiquid
216  (
217  multicomponentThermos.valid().second()
218  ? vfToVif(multicomponentThermos.second().rhoi(specieis().second(), p, T))
219  : vfToVif(thermos().second().rho())
220  );
221  DebugField(WGas);
222  DebugField(rhoGas);
223  DebugField(WLiquid);
224  DebugField(rhoLiquid);
225 
226  // Surface tension
228  (
229  fluid_.sigma
230  (
232  (
233  fluid_.phases()[phaseNames().first()],
234  fluid_.phases()[phaseNames().second()]
235  )
236  )
237  );
238  DebugField(sigma);
239 
240  // Mole fraction of nucleating specie
241  const volScalarField::Internal Xi
242  (
243  thermoGas.Y()[specieis().first()]*WGas/thermoGas.Wi(specieis().first())
244  );
245 
246  // Saturation pressure and concentration
247  const volScalarField::Internal pSat(saturationModel_->pSat(T()));
248  const volScalarField::Internal cSat(pSat/p()*rhoGas/WGas);
249  DebugField(pSat);
250  DebugField(cSat);
251 
252  // Supersaturation of the nucleating specie
253  const volScalarField::Internal S(Xi*p()/pSat);
254  DebugField(S);
255 
256  // Mass, volume and diameter of one molecule in the condensed phase
257  const volScalarField::Internal mMolc(WLiquid/NNA);
258  const volScalarField::Internal vMolc(mMolc/rhoLiquid);
259  const volScalarField::Internal dMolc(cbrt(6/pi*vMolc));
260  DebugField(mMolc);
261  DebugField(vMolc);
262  DebugField(dMolc);
263 
264  // Diameter of nuclei
265  d_ = 4*sigma*vMolc/(k*T()*log(max(S, 1 + small)));
266  DebugField(d_);
267 
268  // ?
269  const volScalarField::Internal deltaPhiStar(pi/3*sigma*sqr(d_));
270  DebugField(deltaPhiStar);
271 
272  // Ratio of nucleus volume to molecular volume
273  const volScalarField::Internal iStar(pi/6*pow3(d_)/vMolc);
274  DebugField(iStar);
275 
276  // ?
277  const volScalarField::Internal betaIStar1
278  (
279  sqrt(6*k*T()/mMolc)*sqrt((iStar + 1)/iStar)*sqr(d_/2 + dMolc/2)
280  );
281  DebugField(betaIStar1);
282 
283  // Number-based nucleation rate; i.e., number of nuclei created per second
284  // per unit volume
286  (
287  betaIStar1*sqr(cSat*NNA)*exp(-deltaPhiStar/(k*T()))
288  *sqrt(sigma/(k*T()))
289  *2*mMolc/(pi*sqr(d_)*rhoLiquid)
290  );
291  DebugField(J);
292 
293  // Mass transfer rate
294  mDotByAlphaGas_ = J*iStar*mMolc;
295  DebugField(mDotByAlphaGas_);
296 
297  #undef DebugField
298 }
299 
300 
302 (
303  const volScalarField& alpha,
304  const volScalarField& rho,
305  fvMatrix<scalar>& eqn
306 ) const
307 {
308  const label i = index(alphaNames(), eqn.psi().name());
309 
310  // !!! Note at present multiphaseEuler cannot linearise w.r.t alphaA in the
311  // continuity equation for alphaB. So we can only create a linearised
312  // source for this model in the gas volume-fraction equation.
313 
314  if (i == 0)
315  {
316  eqn -= fvm::Sp(mDotByAlphaGas_, eqn.psi());
317  }
318  else
319  {
321  }
322 }
323 
324 
326 {
327  if (phaseChange::read(dict))
328  {
329  readCoeffs(coeffs(dict));
330  return true;
331  }
332  else
333  {
334  return false;
335  }
336 }
337 
338 
339 // ************************************************************************* //
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, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< 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
const Type & second() const
Return second.
Definition: PairI.H:121
Class containing a pair of thermo references. Handles down-casting to more specific thermo types by c...
Definition: ThermoRefPair.H:51
const Pair< bool > & valid() const
Access the validity flags.
const ThermoType & first() const
Access the first thermo.
const ThermoType & second() const
Access the second thermo.
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
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.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Use phaseChange's source functions.
Definition: phaseChange.C:551
virtual tmp< DimensionedField< scalar, volMesh > > d() const
Return the diameter of nuclei.
virtual tmp< DimensionedField< scalar, volMesh > > tau() const
Return the nucleation time scale.
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:225
Base class for phase change models.
Definition: phaseChange.H:61
void reReadSpecie(const dictionary &dict) const
Re-read the names of the transferring specie.
Definition: phaseChange.C:129
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
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].
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Class to represent a system of phases.
Definition: phaseSystem.H:74
static autoPtr< saturationPressureModel > New(const dictionary &dict)
Select with dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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))
const dimensionedScalar k
Boltzmann constant.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar NNA
Avagadro number: default SI units: [1/kmol].
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
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)
const dimensionSet dimDensity
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
labelList fv(nPoints)
dictionary dict
volScalarField & p