Kunz.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-2021 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 "Kunz.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace twoPhaseChangeModels
34 {
35  defineTypeNameAndDebug(Kunz, 0);
36  addToRunTimeSelectionTable(cavitationModel, Kunz, dictionary);
37 }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
45 )
46 :
47  cavitationModel(typeName, mixture),
48 
49  UInf_("UInf", dimVelocity, twoPhaseChangeModelCoeffs_),
50  tInf_("tInf", dimTime, twoPhaseChangeModelCoeffs_),
51  Cc_("Cc", dimless, twoPhaseChangeModelCoeffs_),
52  Cv_("Cv", dimless, twoPhaseChangeModelCoeffs_),
53 
54  p0_("0", pSat().dimensions(), 0.0),
55 
56  mcCoeff_(Cc_*mixture_.rho2()/tInf_),
57  mvCoeff_(Cv_*mixture_.rho2()/(0.5*mixture_.rho1()*sqr(UInf_)*tInf_))
58 {
59  correct();
60 }
61 
62 
63 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
64 
67 {
68  const volScalarField& p =
69  mixture_.U().db().lookupObject<volScalarField>("p");
70 
71  const volScalarField limitedAlpha1
72  (
73  min(max(mixture_.alpha1(), scalar(0)), scalar(1))
74  );
75 
77  (
78  mcCoeff_*sqr(limitedAlpha1)
79  *max(p - pSat(), p0_)/max(p - pSat(), 0.01*pSat()),
80 
81  mvCoeff_*min(p - pSat(), p0_)
82  );
83 }
84 
85 
88 {
89  const volScalarField& p =
90  mixture_.U().db().lookupObject<volScalarField>("p");
91 
92  const volScalarField limitedAlpha1
93  (
94  min(max(mixture_.alpha1(), scalar(0)), scalar(1))
95  );
96 
98  (
99  mcCoeff_*sqr(limitedAlpha1)*(1.0 - limitedAlpha1)
100  *pos0(p - pSat())/max(p - pSat(), 0.01*pSat()),
101 
102  (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
103  );
104 }
105 
106 
108 {
110 }
111 
112 
114 {
115  if (cavitationModel::read())
116  {
117  twoPhaseChangeModelCoeffs_ = optionalSubDict(type() + "Coeffs");
118 
119  twoPhaseChangeModelCoeffs_.lookup("UInf") >> UInf_;
120  twoPhaseChangeModelCoeffs_.lookup("tInf") >> tInf_;
121  twoPhaseChangeModelCoeffs_.lookup("Cc") >> Cc_;
122  twoPhaseChangeModelCoeffs_.lookup("Cv") >> Cv_;
123 
124  mcCoeff_ = Cc_*mixture_.rho2()/tInf_;
125  mvCoeff_ = Cv_*mixture_.rho2()/(0.5*mixture_.rho1()*sqr(UInf_)*tInf_);
126 
127  return true;
128  }
129  else
130  {
131  return false;
132  }
133 }
134 
135 
136 // ************************************************************************* //
defineTypeNameAndDebug(cavitationModel, 0)
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
Definition: Kunz.C:87
virtual void correct()
Correct the Kunz phaseChange model.
Definition: Kunz.C:107
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Abstract base class for cavitation models.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Kunz(const immiscibleIncompressibleTwoPhaseMixture &mixture)
Construct for mixture.
Definition: Kunz.C:43
const dimensionSet dimless
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const dimensionSet dimTime
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
Definition: Kunz.C:66
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool read()
Read the transportProperties dictionary and update.
Definition: Kunz.C:113
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual void correct()=0
Correct the phaseChange model.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimVelocity
addToRunTimeSelectionTable(cavitationModel, Kunz, dictionary)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
An immiscible incompressible two-phase mixture transport model.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField & p
Namespace for OpenFOAM.
virtual bool read()
Read the transportProperties dictionary and update.