Kunz.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 "Kunz.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace phaseChangeTwoPhaseMixtures
37 {
38  defineTypeNameAndDebug(Kunz, 0);
39  addToRunTimeSelectionTable(phaseChangeTwoPhaseMixture, Kunz, components);
40 }
41 }
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const volVectorField& U,
48  const surfaceScalarField& phi
49 )
50 :
51  phaseChangeTwoPhaseMixture(typeName, U, phi),
52 
53  UInf_("UInf", dimVelocity, phaseChangeTwoPhaseMixtureCoeffs_),
54  tInf_("tInf", dimTime, phaseChangeTwoPhaseMixtureCoeffs_),
55  Cc_("Cc", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
56  Cv_("Cv", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
57 
58  p0_(pSat().dimensions(), Zero),
59 
60  mcCoeff_(Cc_*rho2()/tInf_),
61  mvCoeff_(Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_))
62 {
63  correct();
64 }
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
71 {
73  volScalarField limitedAlpha1(clamp(alpha1_, zero_one{}));
74 
75  return Pair<tmp<volScalarField>>
76  (
77  mcCoeff_*sqr(limitedAlpha1)
78  *max(p - pSat(), p0_)/max(p - pSat(), 0.01*pSat()),
79 
80  mvCoeff_*min(p - pSat(), p0_)
81  );
82 }
83 
86 {
87  const volScalarField& p = alpha1_.db().lookupObject<volScalarField>("p");
88  volScalarField limitedAlpha1(clamp(alpha1_, zero_one{}));
89 
90  return Pair<tmp<volScalarField>>
91  (
92  mcCoeff_*sqr(limitedAlpha1)*(1.0 - limitedAlpha1)
93  *pos0(p - pSat())/max(p - pSat(), 0.01*pSat()),
94 
95  (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
96  );
97 }
98 
99 
101 {}
102 
103 
105 {
107  {
108  phaseChangeTwoPhaseMixtureCoeffs_ = optionalSubDict(type() + "Coeffs");
109 
110  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("UInf", UInf_);
111  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("tInf", tInf_);
112  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cc", Cc_);
113  phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cv", Cv_);
114 
115  mcCoeff_ = Cc_*rho2()/tInf_;
116  mvCoeff_ = Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_);
117 
118  return true;
119  }
120 
121  return false;
122 }
123 
124 
125 // ************************************************************************* //
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
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 Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & rho1
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
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
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
volScalarField alpha1_
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
const dimensionedScalar & pSat() const
Return const-access to the saturation vapour pressure.
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)
Kunz(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & rho2
U
Definition: pEqn.H:72
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void correct()
Correct the Kunz phaseChange model.
virtual bool read()
Read the transportProperties dictionary and update.
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity