homogeneousLiquidPhaseSeparation.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"
30 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
41  (
42  fvModel,
45  );
46 }
47 }
48 
49 const Foam::NamedEnum
50 <
52  3
53 >
55 {"solid", "liquid", "gas"};
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 void Foam::fv::homogeneousLiquidPhaseSeparation::readCoeffs
61 (
62  const dictionary& dict
63 )
64 {
66 
67  solubilityCurve_.reset
68  (
70  (
71  "solubility",
74  dict
75  ).ptr()
76  );
77 
78  nucleateType_ = nucleateTypeNames_.read(dict.lookup("nucleate"));
79 }
80 
81 
83 Foam::fv::homogeneousLiquidPhaseSeparation::YSat
84 (
86 ) const
87 {
88  return
90  (
91  name() + ":YSat",
92  mesh(),
93  dimless,
94  solubilityCurve_->value(T)
95  );
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
102 (
103  const word& name,
104  const word& modelType,
105  const fvMesh& mesh,
106  const dictionary& dict
107 )
108 :
109  phaseChange(name, modelType, mesh, dict, readSpecie(dict, true)),
110  fluid_
111  (
112  mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)
113  ),
114  d_
115  (
116  IOobject
117  (
118  name + ":d",
119  mesh.time().name(),
120  mesh,
121  IOobject::READ_IF_PRESENT,
122  IOobject::AUTO_WRITE
123  ),
124  mesh,
126  ),
127  mDotByAlphaSolution_
128  (
129  IOobject
130  (
131  name + ":mDotByAlpha",
132  mesh.time().name(),
133  mesh,
134  IOobject::READ_IF_PRESENT,
135  IOobject::AUTO_WRITE
136  ),
137  mesh,
139  ),
140  solubilityCurve_(nullptr),
141  nucleateType_(nucleateType::solid)
142 {
143  readCoeffs(coeffs(dict));
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
151 {
152  return d_;
153 }
154 
155 
158 {
159  const volScalarField::Internal& alphaSolution =
160  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
161 
162  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
163  this->multicomponentThermos(true, false);
164 
165  const multicomponentThermo& thermoSolution = multicomponentThermos.first();
166 
167  const volScalarField& p = this->p();
168  const volScalarField& T = thermoSolution.T();
169 
170  const volScalarField::Internal rhoPrecipitate
171  (
172  multicomponentThermos.valid().second()
173  ? vfToVif(multicomponentThermos.second().rhoi(specieis().second(), p, T))
174  : vfToVif(thermos().second().rho())
175  );
176 
178 
179  return alphaSolution*mDotByAlphaSolution_/(rhoPrecipitate*v);
180 }
181 
182 
185 {
186  const volScalarField::Internal& alphaSolution =
187  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
188 
189  return alphaSolution*mDotByAlphaSolution_;
190 }
191 
192 
195 {
196  static const dimensionedScalar mDotRootVSmall
197  (
199  rootVSmall
200  );
201 
202  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
203  this->multicomponentThermos(true, false);
204 
205  const multicomponentThermo& thermoSolution = multicomponentThermos.first();
206 
207  // Solution density
208  const volScalarField::Internal rhoSolution(vfToVif(thermoSolution.rho()));
209 
210  // Mass fraction of nucleating specie
211  const volScalarField::Internal Yi = thermoSolution.Y()[specieis().first()];
212 
213  return Yi*rhoSolution/max(mDotByAlphaSolution_, mDotRootVSmall);
214 }
215 
216 
218 {
219  #define DebugField(field) \
220  DebugInfo \
221  << name() << ": " \
222  << #field << ' ' << field.dimensions() << " min/avg/max = " \
223  << gMin(field) << '/' << gAverage(field) << '/' << gMax(field) \
224  << nl;
225 
229 
230  const ThermoRefPair<multicomponentThermo> multicomponentThermos =
231  this->multicomponentThermos(true, false);
232 
233  const fluidMulticomponentThermo& thermoSolution =
234  this->fluidMulticomponentThermos(true, false).first();
235 
236  const volScalarField& p = this->p();
237  const volScalarField& T = thermoSolution.T();
238  DebugField(p);
239  DebugField(T);
240 
241  // Phase molecular masses and densities
242  const volScalarField::Internal rhoSolution(vfToVif(thermoSolution.rho()));
243  const volScalarField::Internal muSolution(vfToVif(thermoSolution.mu()));
244  const volScalarField::Internal WPrecipitate
245  (
246  multicomponentThermos.valid().second()
248  (
249  "W",
250  mesh(),
251  multicomponentThermos.second().Wi(specieis().second())
252  )
253  : vfToVif(thermos().second().W())
254  );
255  const volScalarField::Internal rhoPrecipitate
256  (
257  multicomponentThermos.valid().second()
258  ? vfToVif(multicomponentThermos.second().rhoi(specieis().second(), p, T))
259  : vfToVif(thermos().second().rho())
260  );
261  DebugField(rhoSolution);
262  DebugField(WPrecipitate);
263  DebugField(rhoPrecipitate);
264 
265  // Surface tension
267  (
268  fluid_.sigma
269  (
271  (
272  fluid_.phases()[phaseNames().first()],
273  fluid_.phases()[phaseNames().second()]
274  )
275  )
276  );
277  DebugField(sigma);
278 
279  // Mass fraction of nucleating specie
280  const volScalarField::Internal Yi = thermoSolution.Y()[specieis().first()];
281 
282  // Saturation mass fraction and concentration
283  const volScalarField::Internal solubility
284  (
286  (
287  "YSat",
288  mesh(),
289  dimless,
290  solubilityCurve_->value(T)
291  )
292  );
293  const volScalarField::Internal YSat(solubility/(solubility + 1));
294  const volScalarField::Internal cSat(YSat*rhoSolution/WPrecipitate);
295  DebugField(YSat);
296  DebugField(cSat);
297 
298  // Supersaturation of the nucleating specie
299  const volScalarField::Internal S(Yi/YSat);
300  DebugField(S);
301 
302  // Mass and volume of one molecule in the precipitate
303  const volScalarField::Internal mMolc(WPrecipitate/NNA);
304  const volScalarField::Internal vMolc(mMolc/rhoPrecipitate);
305  const volScalarField::Internal dMolc(cbrt(6/pi*vMolc));
306  DebugField(mMolc);
307  DebugField(vMolc);
308  DebugField(dMolc);
309 
310  // Diameter of nuclei
311  d_ = 4*sigma*vMolc/(k*T()*log(max(S, 1 + small)));
312  DebugField(d_);
313 
314  // ?
315  const volScalarField::Internal deltaPhiStar(pi/3*sigma*sqr(d_));
316  DebugField(deltaPhiStar);
317 
318  // Ratio of nucleus volume to molecular volume
319  const volScalarField::Internal iStar(pi/6*pow3(d_)/vMolc);
320  DebugField(iStar);
321 
322  // Pre-exponential factor. Depends on the type of nucleates.
324  switch (nucleateType_)
325  {
326  case nucleateType::solid:
327  case nucleateType::liquid:
328  talpha = k*T()/(3*pi*pow3(dMolc)*muSolution);
329  break;
330 
331  case nucleateType::gas:
332  talpha = sqrt(2*sigma/(pi*mMolc));
333  break;
334  }
335 
336  // Number-based nucleation rate; i.e., number of nuclei created per second
337  // per unit volume
339  (
340  cSat*NNA*talpha*exp(-deltaPhiStar/(k*T()))
341  );
342  DebugField(J);
343 
344  // Mass transfer rate
345  mDotByAlphaSolution_ = J*iStar*mMolc;
346  DebugField(mDotByAlphaSolution_);
347 
348  #undef DebugField
349 }
350 
351 
353 (
354  const volScalarField& alpha,
355  const volScalarField& rho,
356  fvMatrix<scalar>& eqn
357 ) const
358 {
359  const label i = index(alphaNames(), eqn.psi().name());
360 
361  // !!! Note at present multiphaseEuler cannot linearise w.r.t alphaA in the
362  // continuity equation for alphaB. So we can only create a linearised
363  // source for this model in the solution volume-fraction equation.
364 
365  if (i == 0)
366  {
367  eqn -= fvm::Sp(mDotByAlphaSolution_, eqn.psi());
368  }
369  else
370  {
372  }
373 }
374 
375 
377 {
378  if (phaseChange::read(dict))
379  {
380  readCoeffs(coeffs(dict));
381  return true;
382  }
383  else
384  {
385  return false;
386  }
387 }
388 
389 
390 // ************************************************************************* //
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,.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:54
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 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
Base-class for multi-component fluid thermodynamic properties.
virtual const volScalarField & mu() const =0
Dynamic viscosity of mixture [kg/m/s].
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 a solid or liquid phase separating out of a liquid solution.
static const NamedEnum< nucleateType, 3 > nucleateTypeNames_
Names of the moment types.
homogeneousLiquidPhaseSeparation(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
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.
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
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
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimLength
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
const dimensionSet dimTemperature
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)
const unitConversion unitFraction
labelList fv(nPoints)
dictionary dict
volScalarField & p