Last Updated: June 9, 2019
The cavitatingFoam is a transient cavitation solver based on the homogeneous equilibrium model (HEM) from which the compressibility of the liquid/vapour ‘mixture’ is obtained, whose density varies from liquid density to vapor one according to the chosen barotropic equation of state.
Compressibility and speed of sound
\begin{align}
\frac{D \rho}{Dt} = \psi \frac{Dp}{Dt} \tag{1} \label{eq:eos}
\end{align}
\begin{align}
\psi = \rho \beta_s \tag{2} \label{eq:psi}
\end{align}
- \(\rho\): density
- \(p\): pressure
- \(\beta_s\): isentropic compressibility
The definition of the isentropic compressibility can be deformed in the following way:
\begin{align}
\beta_s &\equiv -\frac{1}{V} \left( \frac{\partial V}{\partial p} \right)_s \tag{3} \label{eq:betas} \\
&= -\frac{1}{V} \left(\frac{\partial V}{\partial \rho} \frac{\partial \rho}{\partial p} \right)_s \\
&= \rho \frac{1}{\rho^2} \left(\frac{\partial \rho}{\partial p} \right)_s \\
&= \frac{1}{\rho c^2}
\end{align}
where \(c\) is the speed of sound. We can substitute this expression for \(\beta_s\) in \eqref{eq:psi}, obtaining the following expression:
\begin{align}
\psi = \frac{1}{c^2}. \tag{4} \label{eq:psi-sos}
\end{align}
The compressibility of the mixture equals to the reciprocal of the square of the speed of sound.
We can specify a compressibility function \(\psi\) that defines the coupling between the density and pressure in the thermophysicalProperties file.
barotropicCompressibilityModel
The following three models for \(\psi\) are available in OpenFOAM v1812.
- linear
- Chung
- Wallis
| 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 | /*--------------------------------*- C++ -*----------------------------------*\ | =========                 |                                                 | | \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           | |  \\    /   O peration     | Version:  v1812                                 | |   \\  /    A nd           | Web:      www.OpenFOAM.com                      | |    \\/     M anipulation  |                                                 | \*---------------------------------------------------------------------------*/ FoamFile {     version     2.0;     format      ascii;     class       dictionary;     location    "constant";     object      thermodynamicProperties; } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // barotropicCompressibilityModel linear; psiv            2.5e-06; rholSat         830; psil            5e-07; pSat            4500; rhoMin          0.001; // ************************************************************************* // | 
- psiv: compressibility \(\psi\) of the vapor phase
- psil: compressibility \(\psi\) of the liquid phase
- pSat: saturation pressure
- rhovSat: density of vapor at saturation point
- rholSat: density of liquid at saturation point
The expressions of the speed of sound \(c\) in the liquid/vapour mixture are different according to the selected function.
linear
\begin{align}
\psi = \gamma \psi_v + \left(1 – \gamma \right)\psi_l \tag{5} \label{eq:linear}
\end{align}
| 75 76 77 78 | void Foam::compressibilityModels::linear::correct() {     psi_ = gamma_*psiv_ + (scalar(1) - gamma_)*psil_; } | 
Chung
\begin{align}
s_{fa} = \sqrt{ \frac{\frac{\rho_{v, Sat}}{\psi_v}}{\left( 1- \gamma \right)\frac{\rho_{v, Sat}}{\psi_v} + \gamma \frac{\rho_{l, Sat}}{\psi_l}} }  \tag{6} \label{eq:chung-sfa}
\end{align}
\begin{align}
\psi = \left( \left( \frac{1 – \gamma}{\sqrt{\psi_v}} + \gamma \frac{s_{fa}}{\sqrt{\psi_l}} \right)\frac{\sqrt{\psi_v \psi_l}}{s_{fa}} \right)^2 \tag{7} \label{eq:chung}
\end{align}
| 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 | void Foam::compressibilityModels::Chung::correct() {     volScalarField sfa     (         sqrt         (             (rhovSat_/psiv_)            /((scalar(1) - gamma_)*rhovSat_/psiv_ + gamma_*rholSat_/psil_)         )     );     psi_ = sqr     (         ((scalar(1) - gamma_)/sqrt(psiv_) + gamma_*sfa/sqrt(psil_))        *sqrt(psiv_*psil_)/sfa     ); } | 
Wallis
\begin{align}
\psi = \left( \gamma \rho_{v, Sat} + \left( 1 – \gamma \right)\rho_{l, Sat} \right)
\left( \gamma \frac{\psi_v}{\rho_{v, Sat}} + \left( 1 – \gamma \right)\frac{\psi_l}{\rho_{l, Sat}} \right) \tag{8} \label{eq:wallis}
\end{align}
\begin{align}
\frac{\psi}{\gamma \rho_{v, Sat} + \left( 1 – \gamma \right)\rho_{l, Sat}} = \gamma \frac{\psi_v}{\rho_{v, Sat}} + \left( 1 – \gamma \right)\frac{\psi_l}{\rho_{l, Sat}} \tag{8′} \label{eq:wallis2}
\end{align}
| 86 87 88 89 90 91 | void Foam::compressibilityModels::Wallis::correct() {     psi_ =         (gamma_*rhovSat_ + (scalar(1) - gamma_)*rholSat_)        *(gamma_*psiv_/rhovSat_ + (scalar(1) - gamma_)*psil_/rholSat_); } | 
References
[1] Christopher Earls Brennen, Fundamentals of Multiphase Flow. Cambridge University Press (2005)
 
	