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 (
44  const compressibleTwoPhaseMixture& mixture
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  correct();
57 }
58 
59 
60 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
61 
64 {
65  const volScalarField::Internal& p =
67 
68  const volScalarField::Internal mcCoeff_(Cc_*rho2()/tInf_);
69  const volScalarField::Internal mvCoeff_
70  (
71  Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_)
72  );
73 
74  const volScalarField::Internal limitedAlpha1
75  (
76  min(max(mixture_.alpha1()(), scalar(0)), scalar(1))
77  );
78 
79  return Pair<tmp<volScalarField::Internal>>
80  (
81  mcCoeff_*sqr(limitedAlpha1)
82  *max(p - pSat(), p0_)/max(p - pSat(), 0.01*pSat()),
83 
84  mvCoeff_*min(p - pSat(), p0_)
85  );
86 }
87 
88 
91 {
92  const volScalarField::Internal& p =
94 
95  const volScalarField::Internal mcCoeff_(Cc_*rho2()/tInf_);
96  const volScalarField::Internal mvCoeff_
97  (
98  Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_)
99  );
100 
101  const volScalarField::Internal limitedAlpha1
102  (
103  min(max(mixture_.alpha1()(), scalar(0)), scalar(1))
104  );
105 
106  return Pair<tmp<volScalarField::Internal>>
107  (
108  mcCoeff_*sqr(limitedAlpha1)*(1.0 - limitedAlpha1)
109  *pos0(p - pSat())/max(p - pSat(), 0.01*pSat()),
110 
111  (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
112  );
113 }
114 
115 
117 {
119 }
120 
121 
123 {
124  if (cavitationModel::read())
125  {
127 
128  twoPhaseChangeModelCoeffs_.lookup("UInf") >> UInf_;
129  twoPhaseChangeModelCoeffs_.lookup("tInf") >> tInf_;
130  twoPhaseChangeModelCoeffs_.lookup("Cc") >> Cc_;
131  twoPhaseChangeModelCoeffs_.lookup("Cv") >> Cv_;
132 
133  return true;
134  }
135  else
136  {
137  return false;
138  }
139 }
140 
141 
142 // ************************************************************************* //
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
dictionary twoPhaseChangeModelCoeffs_
Model coefficient dictionary.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Kunz(const immiscibleIncompressibleTwoPhaseMixture &mixture)
Construct for mixture.
Definition: Kunz.C:43
const dimensionSet dimless
const dimensionedScalar & pSat() const
Return const-access to the saturation vapour pressure.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
Definition: Kunz.C:66
const volScalarField::Internal & rho1() const
Return the internal field of the density of phase 1.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool read()
Read the phaseProperties dictionary and update.
Definition: Kunz.C:113
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1060
virtual void correct()=0
Correct the phaseChange model.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimVelocity
addToRunTimeSelectionTable(cavitationModel, Kunz, dictionary)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
const volScalarField::Internal & rho2() const
Return the internal field of the density of phase 2.
const immiscibleIncompressibleTwoPhaseMixture & mixture_
Reference to the two-phase mixture.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
virtual bool read()
Read the phaseProperties dictionary and update.