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-2025 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 compressible
34 {
35 namespace cavitationModels
36 {
39 }
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const dictionary& dict,
49  const compressibleTwoPhases& phases,
50  const label liquidIndex
51 )
52 :
53  cavitationModel(dict, phases, liquidIndex),
54 
55  UInf_("UInf", dimVelocity, dict),
56  tInf_("tInf", dimTime, dict),
57  Cv_("Cv", dimless, dict),
58  Cc_("Cc", dimless, dict),
59 
60  p0_("0", dimPressure, 0)
61 {
62  correct();
63 }
64 
65 
66 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
67 
69 Foam::compressible::cavitationModels::Kunz::mvCoeff() const
70 {
71  return Cv_*rhov()/(0.5*rhol()*sqr(UInf_)*tInf_);
72 }
73 
74 
76 Foam::compressible::cavitationModels::Kunz::mcCoeff() const
77 {
78  return Cc_*rhov()/tInf_;
79 }
80 
81 
84 {
86  phases_.mesh().lookupObject<volScalarField>("p");
87 
88  const volScalarField::Internal alphal
89  (
90  min(max(this->alphal(), scalar(0)), scalar(1))
91  );
92 
93  const volScalarField::Internal pSatv(this->pSatv());
94  const volScalarField::Internal pSatl(this->pSatl());
95 
97  (
98  mcCoeff()*sqr(alphal)
99  *max(p - pSatv, p0_)/max(p - pSatv, 0.01*pSatv),
100  - mvCoeff()*min(p - pSatl, p0_)
101  );
102 }
103 
104 
107 {
108  const volScalarField::Internal& p =
109  phases_.mesh().lookupObject<volScalarField>("p");
110 
111  const volScalarField::Internal alphav
112  (
113  min(max(this->alphav(), scalar(0)), scalar(1))
114  );
115 
116  const volScalarField::Internal alphal
117  (
118  min(max(this->alphal(), scalar(0)), scalar(1))
119  );
120 
121  const volScalarField::Internal pSatv(this->pSatv());
122  const volScalarField::Internal pSatl(this->pSatl());
123 
125  (
126  mcCoeff()*alphav*sqr(alphal)*pos0(p - pSatv)
127  /max(p - pSatv, 0.01*pSatv),
128  - mvCoeff()*alphal*neg(p - pSatl)
129  );
130 }
131 
132 
134 {}
135 
136 
138 (
139  const dictionary& dict
140 )
141 {
143  {
144  dict.lookup("UInf") >> UInf_;
145  dict.lookup("tInf") >> tInf_;
146  dict.lookup("Cv") >> Cv_;
147  dict.lookup("Cc") >> Cc_;
148 
149  return true;
150  }
151  else
152  {
153  return false;
154  }
155 }
156 
157 
158 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
virtual bool read(const dictionary &dict)=0
Read the dictionary and update.
virtual Pair< tmp< volScalarField::Internal > > mDotcvAlphal() const
Return the mass condensation and vaporisation rates as a.
Definition: Kunz.C:83
virtual void correct()
Correct the Kunz phaseChange model.
Definition: Kunz.C:133
Kunz(const dictionary &dict, const compressibleTwoPhases &phases, const label liquidIndex)
Construct for phases.
Definition: Kunz.C:47
virtual bool read(const dictionary &dict)
Read the dictionary and update.
Definition: Kunz.C:138
virtual Pair< tmp< volScalarField::Internal > > mDotcvP() const
Return the mass condensation and vaporisation rates as coefficients.
Definition: Kunz.C:106
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
addToRunTimeSelectionTable(cavitationModel, Kunz, dictionary)
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
dictionary dict
volScalarField & p