Chung.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-2018 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 "Chung.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  namespace compressibilityModels
34  {
35  defineTypeNameAndDebug(Chung, 0);
37  (
38  barotropicCompressibilityModel,
39  Chung,
40  dictionary
41  );
42  }
43 }
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const dictionary& compressibilityProperties,
50  const volScalarField& gamma,
51  const word& psiName
52 )
53 :
54  barotropicCompressibilityModel(compressibilityProperties, gamma, psiName),
55  psiv_
56  (
57  "psiv",
59  compressibilityProperties_.lookup("psiv")
60  ),
61  psil_
62  (
63  "psil",
65  compressibilityProperties_.lookup("psil")
66  ),
67  rhovSat_
68  (
69  "rhovSat",
70  dimDensity,
71  compressibilityProperties_.lookup("rhovSat")
72  ),
73  rholSat_
74  (
75  "rholSat",
76  dimDensity,
77  compressibilityProperties_.lookup("rholSat")
78  )
79 {
80  correct();
81 }
82 
83 
84 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
85 
87 {
88  volScalarField sfa
89  (
90  sqrt
91  (
92  (rhovSat_/psiv_)
93  /((scalar(1) - gamma_)*rhovSat_/psiv_ + gamma_*rholSat_/psil_)
94  )
95  );
96 
97  psi_ = sqr
98  (
99  ((scalar(1) - gamma_)/sqrt(psiv_) + gamma_*sfa/sqrt(psil_))
100  *sqrt(psiv_*psil_)/sfa
101  );
102 }
103 
104 
106 (
107  const dictionary& compressibilityProperties
108 )
109 {
110  barotropicCompressibilityModel::read(compressibilityProperties);
111 
112  compressibilityProperties_.lookup("psiv") >> psiv_;
113  compressibilityProperties_.lookup("psil") >> psil_;
114  compressibilityProperties_.lookup("rhovSat") >> rhovSat_;
115  compressibilityProperties_.lookup("rholSat") >> rholSat_;
116 
117  return true;
118 }
119 
120 
121 // ************************************************************************* //
bool read(const dictionary &compressibilityProperties)
Read transportProperties dictionary.
Definition: Chung.C:106
virtual bool read(const dictionary &compressibilityProperties)=0
Read compressibilityProperties dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void correct()
Correct the Chung compressibility.
Definition: Chung.C:86
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/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:72
addToRunTimeSelectionTable(barotropicCompressibilityModel, Chung, dictionary)
const dimensionSet dimDensity
Abstract class for barotropic compressibility models.
const dimensionSet dimCompressibility
Namespace for OpenFOAM.
Chung(const dictionary &compressibilityProperties, const volScalarField &gamma, const word &psiName="psi")
Construct from components.
Definition: Chung.C:48