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 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::homogeneousLiquidPhaseSeparation::readCoeffs()
52 {
53  solubilityCurve_.reset
54  (
56  (
57  "solubility",
60  coeffs()
61  ).ptr()
62  );
63 }
64 
65 
67 Foam::fv::homogeneousLiquidPhaseSeparation::YSat
68 (
70 ) const
71 {
72  return
74  (
75  typedName(IOobject::groupName("YSat", interface_.name())),
76  mesh(),
77  dimless,
78  solubilityCurve_->value(T)
79  );
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 (
87  const word& name,
88  const word& modelType,
89  const fvMesh& mesh,
90  const dictionary& dict
91 )
92 :
94  (
95  name,
96  modelType,
97  mesh,
98  dict,
99  {true, false},
100  {true, false}
101  ),
102  fluid_
103  (
104  mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)
105  ),
106  interface_
107  (
108  fluid_.phases()[phaseNames().first()],
109  fluid_.phases()[phaseNames().second()]
110  ),
111  d_
112  (
113  IOobject
114  (
115  typedName(IOobject::groupName("d", interface_.name())),
116  mesh.time().name(),
117  mesh,
120  ),
121  mesh,
123  ),
124  mDotByAlphaSolution_
125  (
126  IOobject
127  (
128  typedName(IOobject::groupName("mDotByAlpha", interface_.name())),
129  mesh.time().name(),
130  mesh,
133  ),
134  mesh,
136  ),
137  solubilityCurve_(nullptr)
138 {
139  readCoeffs();
140 }
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
147 {
148  return d_;
149 }
150 
151 
154 {
155  const volScalarField::Internal& alphaSolution =
156  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
157 
158  const multicomponentThermo& thermoSolution = specieThermos().first();
159 
160  const volScalarField& p = this->p();
161  const volScalarField& T = thermoSolution.T();
162 
163  const volScalarField::Internal rhoPrecipitate
164  (
165  specieThermos().valid().second()
166  ? vfToVif(specieThermos().second().rhoi(specieis().second(), p, T))
167  : vfToVif(thermos().second().rho())
168  );
169 
171 
172  return alphaSolution*mDotByAlphaSolution_/(rhoPrecipitate*v);
173 }
174 
175 
178 {
179  const volScalarField::Internal& alphaSolution =
180  mesh().lookupObject<volScalarField::Internal>(alphaNames().first());
181 
182  return alphaSolution*mDotByAlphaSolution_;
183 }
184 
185 
187 {
188  #define DebugField(field) \
189  DebugInfo \
190  << name() << ": " \
191  << #field << ' ' << field.dimensions() << " min/avg/max = " \
192  << gMin(field) << '/' << gAverage(field) << '/' << gMax(field) \
193  << nl;
194 
198 
199  const fluidMulticomponentThermo& thermoSolution =
200  refCast<const fluidMulticomponentThermo>(specieThermos().first());
201 
202  const volScalarField& p = this->p();
203  const volScalarField& T = thermoSolution.T();
204  DebugField(p);
205  DebugField(T);
206 
207  // Phase molecular masses and densities
208  const volScalarField::Internal rhoSolution(vfToVif(thermoSolution.rho()));
209  const volScalarField::Internal muSolution(vfToVif(thermoSolution.mu()));
210  const volScalarField::Internal WPrecipitate
211  (
212  specieThermos().valid().second()
214  (
215  "W",
216  mesh(),
217  specieThermos().second().Wi(specieis().second())
218  )
219  : vfToVif(thermos().second().W())
220  );
221  const volScalarField::Internal rhoPrecipitate
222  (
223  specieThermos().valid().second()
224  ? vfToVif(specieThermos().second().rhoi(specieis().second(), p, T))
225  : vfToVif(thermos().second().rho())
226  );
227  DebugField(rhoSolution);
228  DebugField(WPrecipitate);
229  DebugField(rhoPrecipitate);
230 
231  // Surface tension
232  const volScalarField::Internal sigma(interface_.sigma());
233  DebugField(sigma);
234 
235  // Mass fraction of nucleating specie
236  const volScalarField::Internal Yi = thermoSolution.Y()[specieis().first()];
237 
238  // Saturation mass fraction and concentration
239  const volScalarField::Internal solubility
240  (
242  (
243  "YSat",
244  mesh(),
245  dimless,
246  solubilityCurve_->value(T)
247  )
248  );
249  const volScalarField::Internal YSat(solubility/(solubility + 1));
250  const volScalarField::Internal cSat(YSat*rhoSolution/WPrecipitate);
251  DebugField(YSat);
252  DebugField(cSat);
253 
254  // Supersaturation of the nucleating specie
255  const volScalarField::Internal S(Yi/YSat);
256  DebugField(S);
257 
258  // Mass and volume of one molecule in the precipitate
259  const volScalarField::Internal mMolc(WPrecipitate/NA);
260  const volScalarField::Internal vMolc(mMolc/rhoPrecipitate);
261  const volScalarField::Internal dMolc(cbrt(6/pi*vMolc));
262  DebugField(mMolc);
263  DebugField(vMolc);
264  DebugField(dMolc);
265 
266  // Diameter of nuclei
267  d_ = 4*sigma*vMolc/(k*T()*log(max(S, 1 + small)));
268  DebugField(d_);
269 
270  // ?
271  const volScalarField::Internal deltaPhiStar(pi/3*sigma*sqr(d_));
272  DebugField(deltaPhiStar);
273 
274  // Ratio of nucleus volume to molecular volume
275  const volScalarField::Internal iStar(pi/6*pow3(d_)/vMolc);
276  DebugField(iStar);
277 
278  // Number-based nucleation rate; i.e., number of nuclei created per second
279  // per unit volume
280  const volScalarField::Internal intermediate(-deltaPhiStar/(k*T()));
281  DebugField(intermediate);
283  (
284  cSat*NA*exp(-deltaPhiStar/(k*T()))*k*T()/(3*pi*pow3(dMolc)*muSolution)
285  );
286  DebugField(J);
287 
288  // Mass transfer rate
289  mDotByAlphaSolution_ = J*iStar*mMolc;
290  DebugField(mDotByAlphaSolution_);
291 
292  #undef DebugField
293 }
294 
295 
297 (
298  const volScalarField& alpha,
299  const volScalarField& rho,
300  fvMatrix<scalar>& eqn
301 ) const
302 {
303  const label i = index(alphaNames(), eqn.psi().name());
304 
305  if (i != -1)
306  {
307  if (i == 0)
308  {
309  eqn -= fvm::Sp(mDotByAlphaSolution_, eqn.psi());
310  }
311  else
312  {
313  eqn +=
314  mDot()
315  - correction(fvm::Sp(mDotByAlphaSolution_, eqn.psi()));
316  }
317  }
318  else
319  {
321  }
322 }
323 
324 
326 (
327  const volScalarField& alpha,
328  const volScalarField& rho,
329  const volScalarField& Yi,
330  fvMatrix<scalar>& eqn
331 ) const
332 {
333  const label i = index(alphaNames(), eqn.psi().name());
334 
335  if (i == 0 && specieis().first() != -1 && Yi.member() == specie())
336  {
337  eqn -= fvm::Sp(alpha*mDotByAlphaSolution_, Yi);
338  }
339  else
340  {
342  }
343 }
344 
345 
347 {
349  {
350  readCoeffs();
351  return true;
352  }
353  else
354  {
355  return false;
356  }
357 }
358 
359 
360 // ************************************************************************* //
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,.
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 > 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
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 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
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: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 a solid or liquid phase separating out of a liquid solution.
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.
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.
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
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 dimless
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
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)
const unitConversion unitFraction
labelList fv(nPoints)
dictionary dict
volScalarField & p