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 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 )
57 :
58  cavitationModel(dict, phases),
59 
60  Ca_("Ca", dimless/dimLength, dict),
61  Cv_("Cv", dimless, dict),
62  Cc_("Cc", dimless, dict),
63  alphaNuc_("alphaNuc", dimless, dict),
64 
65  p0_("0", dimPressure, 0)
66 {
67  correct();
68 }
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 Foam::compressible::cavitationModels::Saito::fT
75 (
76  const rhoFluidThermo& thermo
77 ) const
78 {
79  return sqrt(2*pi*(RR/thermo.W()()())*thermo.T()());
80 }
81 
82 
85 {
86  const volScalarField::Internal& p = thermol().p();
87 
88  const volScalarField::Internal alphav
89  (
90  min(max(this->alphav(), scalar(0)), scalar(1))
91  );
92 
93  const volScalarField::Internal alphal
94  (
95  min(max(this->alphal(), scalar(0)), scalar(1))
96  );
97 
98  const volScalarField::Internal alphavNuc(max(alphav, alphaNuc_));
99 
100  const volScalarField::Internal A(Ca_*alphal*alphavNuc);
101 
102  const volScalarField::Internal mvCoeff(Cv_*A*rhol()/(rhov()*fT(thermol())));
103  const volScalarField::Internal mcCoeff(Cc_*A/fT(thermol()));
104 
106  (
107  mcCoeff*alphal*max(p - pSatv(), p0_),
108  -mvCoeff*alphavNuc*min(p - pSatl(), p0_)
109  );
110 }
111 
112 
115 {
116  const volScalarField::Internal& p = thermol().p();
117 
118  const volScalarField::Internal alphav
119  (
120  min(max(this->alphav(), scalar(0)), scalar(1))
121  );
122 
123  const volScalarField::Internal alphal
124  (
125  min(max(this->alphal(), scalar(0)), scalar(1))
126  );
127 
128  const volScalarField::Internal alphavNuc(max(alphav, alphaNuc_));
129 
130  const volScalarField::Internal A(Ca_*alphal*alphavNuc);
131 
132  const volScalarField::Internal mvCoeff(Cv_*A*rhol()/(rhov()*fT(thermol())));
133  const volScalarField::Internal mcCoeff(Cc_*A/fT(thermol()));
134 
136  (
137  mcCoeff*alphal*alphav*pos0(p - pSatv()),
138  -mvCoeff*alphal*alphavNuc*neg(p - pSatl())
139  );
140 }
141 
142 
144 {}
145 
146 
148 (
149  const dictionary& dict
150 )
151 {
153  {
154  dict.lookup("Ca") >> Ca_;
155  dict.lookup("Cv") >> Cv_;
156  dict.lookup("Cc") >> Cc_;
157  dict.lookup("alphaNuc") >> alphaNuc_;
158 
159  return true;
160  }
161  else
162  {
163  return false;
164  }
165 }
166 
167 
168 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
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 <T> with first() and second() elements.
Definition: Pair.H:65
virtual bool read(const dictionary &dict)=0
Read the dictionary and update.
Saito(const dictionary &dict, const compressibleTwoPhases &phases)
Construct for phases.
Definition: Saito.C:53
virtual Pair< tmp< volScalarField::Internal > > mDotcvAlphal() const
Return the mass condensation and vaporisation rates as a.
Definition: Saito.C:84
virtual void correct()
Correct the Saito phaseChange model.
Definition: Saito.C:143
virtual bool read(const dictionary &dict)
Read the dictionary and update.
Definition: Saito.C:148
virtual Pair< tmp< volScalarField::Internal > > mDotcvP() const
Return the mass condensation and vaporisation rates as coefficients.
Definition: Saito.C:114
A list of keyword definitions, which are a keyword followed by any number of values (e....
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].
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimPressure
const dimensionSet dimless
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31