Saito.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) 2023-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 "Saito.H"
28 
29 #include "mathematicalConstants.H"
31 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace compressible
40 {
41 namespace cavitationModels
42 {
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict,
55  const compressibleTwoPhases& phases,
56  const label liquidIndex
57 )
58 :
59  cavitationModel(dict, phases, liquidIndex),
60 
61  Ca_("Ca", dimless/dimLength, dict),
62  Cv_("Cv", dimless, dict),
63  Cc_("Cc", dimless, dict),
64  alphaNuc_("alphaNuc", dimless, dict),
65 
66  p0_("0", dimPressure, 0)
67 {
68  correct();
69 }
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 Foam::compressible::cavitationModels::Saito::fT
76 (
77  const rhoFluidThermo& thermo
78 ) const
79 {
80  return sqrt(2*pi*(RR/thermo.W()()())*thermo.T()());
81 }
82 
83 
86 {
87  const volScalarField::Internal& p = thermol().p();
88 
89  const volScalarField::Internal alphav
90  (
91  min(max(this->alphav(), scalar(0)), scalar(1))
92  );
93 
94  const volScalarField::Internal alphal
95  (
96  min(max(this->alphal(), scalar(0)), scalar(1))
97  );
98 
99  const volScalarField::Internal alphavNuc(max(alphav, alphaNuc_));
100 
101  const volScalarField::Internal A(Ca_*alphal*alphavNuc);
102 
103  const volScalarField::Internal mvCoeff(Cv_*A*rhol()/(rhov()*fT(thermol())));
104  const volScalarField::Internal mcCoeff(Cc_*A/fT(thermol()));
105 
107  (
108  mcCoeff*alphal*max(p - pSatv(), p0_),
109  -mvCoeff*alphavNuc*min(p - pSatl(), p0_)
110  );
111 }
112 
113 
116 {
117  const volScalarField::Internal& p = thermol().p();
118 
119  const volScalarField::Internal alphav
120  (
121  min(max(this->alphav(), scalar(0)), scalar(1))
122  );
123 
124  const volScalarField::Internal alphal
125  (
126  min(max(this->alphal(), scalar(0)), scalar(1))
127  );
128 
129  const volScalarField::Internal alphavNuc(max(alphav, alphaNuc_));
130 
131  const volScalarField::Internal A(Ca_*alphal*alphavNuc);
132 
133  const volScalarField::Internal mvCoeff(Cv_*A*rhol()/(rhov()*fT(thermol())));
134  const volScalarField::Internal mcCoeff(Cc_*A/fT(thermol()));
135 
137  (
138  mcCoeff*alphal*alphav*pos0(p - pSatv()),
139  -mvCoeff*alphal*alphavNuc*neg(p - pSatl())
140  );
141 }
142 
143 
145 {}
146 
147 
149 (
150  const dictionary& dict
151 )
152 {
154  {
155  dict.lookup("Ca") >> Ca_;
156  dict.lookup("Cv") >> Cv_;
157  dict.lookup("Cc") >> Cc_;
158  dict.lookup("alphaNuc") >> alphaNuc_;
159 
160  return true;
161  }
162  else
163  {
164  return false;
165  }
166 }
167 
168 
169 // ************************************************************************* //
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...
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: Saito.C:85
virtual void correct()
Correct the Saito phaseChange model.
Definition: Saito.C:144
virtual bool read(const dictionary &dict)
Read the dictionary and update.
Definition: Saito.C:149
Saito(const dictionary &dict, const compressibleTwoPhases &phases, const label liquidIndex)
Construct for phases.
Definition: Saito.C:53
virtual Pair< tmp< volScalarField::Internal > > mDotcvP() const
Return the mass condensation and vaporisation rates as coefficients.
Definition: Saito.C:115
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base-class for fluid thermodynamic properties based on density.
A class for managing temporary objects.
Definition: tmp.H:55
addToRunTimeSelectionTable(cavitationModel, Kunz, dictionary)
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
static const coefficient A("A", dimPressure, 611.21)
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
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31