SchnerrSauer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "SchnerrSauer.H"
30 #include "mathematicalConstants.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace phaseChangeTwoPhaseMixtures
38 {
39  defineTypeNameAndDebug(SchnerrSauer, 0);
41  (
42  phaseChangeTwoPhaseMixture,
43  SchnerrSauer,
44  components
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const volVectorField& U,
55  const surfaceScalarField& phi
56 )
57 :
58  phaseChangeTwoPhaseMixture(typeName, U, phi),
59 
60  n_("n", dimless/dimVolume, phaseChangeTwoPhaseMixtureCoeffs_),
61  dNuc_("dNuc", dimLength, phaseChangeTwoPhaseMixtureCoeffs_),
62  Cc_("Cc", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
63  Cv_("Cv", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
64 
65  p0_(pSat().dimensions(), Zero)
66 {
67  correct();
68 }
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::rRb
75 (
76  const volScalarField& limitedAlpha1
77 ) const
78 {
79  return pow
80  (
82  *limitedAlpha1/(1.0 + alphaNuc() - limitedAlpha1),
83  1.0/3.0
84  );
85 }
86 
87 
89 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::alphaNuc() const
90 {
92  return Vnuc/(1 + Vnuc);
93 }
94 
95 
97 Foam::phaseChangeTwoPhaseMixtures::SchnerrSauer::pCoeff
98 (
99  const volScalarField& p
100 ) const
101 {
102  volScalarField limitedAlpha1(clamp(alpha1_, zero_one{}));
104  (
105  limitedAlpha1*rho1() + (scalar(1) - limitedAlpha1)*rho2()
106  );
107 
108  return
109  (3*rho1()*rho2())*sqrt(2/(3*rho1()))
110  *rRb(limitedAlpha1)/(rho*sqrt(mag(p - pSat()) + 0.01*pSat()));
111 }
112 
113 
116 {
117  const volScalarField& p = alpha1_.db().lookupObject<volScalarField>("p");
118  volScalarField pCoeff(this->pCoeff(p));
119 
120  volScalarField limitedAlpha1(clamp(alpha1_, zero_one{}));
121 
122  return Pair<tmp<volScalarField>>
123  (
124  Cc_*limitedAlpha1*pCoeff*max(p - pSat(), p0_),
125 
126  Cv_*(1.0 + alphaNuc() - limitedAlpha1)*pCoeff*min(p - pSat(), p0_)
127  );
128 }
129 
130 
133 {
134  const volScalarField& p = alpha1_.db().lookupObject<volScalarField>("p");
135  volScalarField pCoeff(this->pCoeff(p));
136 
137  volScalarField limitedAlpha1(clamp(alpha1_, zero_one{}));
138  volScalarField apCoeff(limitedAlpha1*pCoeff);
139 
140  return Pair<tmp<volScalarField>>
141  (
142  Cc_*(1.0 - limitedAlpha1)*pos0(p - pSat())*apCoeff,
143 
144  (-Cv_)*(1.0 + alphaNuc() - limitedAlpha1)*neg(p - pSat())*apCoeff
145  );
146 }
147 
148 
150 {}
151 
152 
154 {
156  {
157  phaseChangeTwoPhaseMixtureCoeffs_ = optionalSubDict(type() + "Coeffs");
158 
159  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("n", n_);
160  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("dNuc", dNuc_);
161  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cc", Cc_);
162  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cv", Cv_);
163 
164  return true;
165  }
166 
167  return false;
168 }
169 
170 
171 // ************************************************************************* //
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
virtual bool read()=0
Read the transportProperties dictionary and update.
virtual void correct()
Correct the SchnerrSauer phaseChange model.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
volScalarField & rho1
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.
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
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:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:425
constexpr scalar pi(M_PI)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
dimensionedScalar pos0(const dimensionedScalar &ds)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & rho2
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
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.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SchnerrSauer(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133