SchnerrSauer.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) 2011-2018 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 
26 #include "SchnerrSauer.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace phaseChangeTwoPhaseMixtures
35 {
36  defineTypeNameAndDebug(SchnerrSauer, 0);
38  (
39  phaseChangeTwoPhaseMixture,
40  SchnerrSauer,
41  components
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const volVectorField& U,
52  const surfaceScalarField& phi
53 )
54 :
55  phaseChangeTwoPhaseMixture(typeName, U, phi),
56 
57  n_("n", dimless/dimVolume, phaseChangeTwoPhaseMixtureCoeffs_),
58  dNuc_("dNuc", dimLength, phaseChangeTwoPhaseMixtureCoeffs_),
59  Cc_("Cc", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
60  Cv_("Cv", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
61 
62  p0_("0", pSat().dimensions(), 0.0)
63 {
64  correct();
65 }
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
71 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::rRb
72 (
73  const volScalarField& limitedAlpha1
74 ) const
75 {
76  return pow
77  (
79  *limitedAlpha1/(1.0 + alphaNuc() - limitedAlpha1),
80  1.0/3.0
81  );
82 }
83 
84 
86 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::alphaNuc() const
87 {
89  return Vnuc/(1 + Vnuc);
90 }
91 
92 
94 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::pCoeff
95 (
96  const volScalarField& p
97 ) const
98 {
99  volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
101  (
102  limitedAlpha1*rho1() + (scalar(1) - limitedAlpha1)*rho2()
103  );
104 
105  return
106  (3*rho1()*rho2())*sqrt(2/(3*rho1()))
107  *rRb(limitedAlpha1)/(rho*sqrt(mag(p - pSat()) + 0.01*pSat()));
108 }
109 
110 
113 {
115  volScalarField pCoeff(this->pCoeff(p));
116 
117  volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
118 
119  return Pair<tmp<volScalarField>>
120  (
121  Cc_*limitedAlpha1*pCoeff*max(p - pSat(), p0_),
122 
123  Cv_*(1.0 + alphaNuc() - limitedAlpha1)*pCoeff*min(p - pSat(), p0_)
124  );
125 }
126 
127 
130 {
132  volScalarField pCoeff(this->pCoeff(p));
133 
134  volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
135  volScalarField apCoeff(limitedAlpha1*pCoeff);
136 
137  return Pair<tmp<volScalarField>>
138  (
139  Cc_*(1.0 - limitedAlpha1)*pos0(p - pSat())*apCoeff,
140 
141  (-Cv_)*(1.0 + alphaNuc() - limitedAlpha1)*neg(p - pSat())*apCoeff
142  );
143 }
144 
145 
147 {}
148 
149 
151 {
153  {
155 
157  phaseChangeTwoPhaseMixtureCoeffs_.lookup("dNuc") >> dNuc_;
160 
161  return true;
162  }
163  else
164  {
165  return false;
166  }
167 }
168 
169 
170 // ************************************************************************* //
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
virtual bool read()=0
Read the transportProperties dictionary and update.
virtual void correct()
Correct the SchnerrSauer phaseChange model.
surfaceScalarField & phi
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void correct()=0
Correct the phaseChange model.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
Generic dimensioned Type class.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual bool read()
Read the transportProperties dictionary and update.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
volScalarField alpha1_
const dimensionedScalar & pSat() const
Return const-access to the saturation vapour pressure.
dimensionedScalar pos0(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
SchnerrSauer(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const dimensionedScalar & rho1() const
Return const-access to phase1 density.