Saito cavitation model. More...
Public Member Functions | |
TypeName ("Saito") | |
Runtime type information. More... | |
Saito (const dictionary &dict, const compressibleTwoPhases &phases) | |
Construct for phases. More... | |
virtual | ~Saito () |
Destructor. More... | |
virtual Pair< tmp< volScalarField::Internal > > | mDotcvAlphal () const |
Return the mass condensation and vaporisation rates as a. More... | |
virtual Pair< tmp< volScalarField::Internal > > | mDotcvP () const |
Return the mass condensation and vaporisation rates as coefficients. More... | |
virtual void | correct () |
Correct the Saito phaseChange model. More... | |
virtual bool | read (const dictionary &dict) |
Read the dictionary and update. More... | |
Public Member Functions inherited from cavitationModel | |
TypeName ("cavitation") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, cavitationModel, dictionary,(const dictionary &dict, const compressibleTwoPhases &phases),(dict, phases)) | |
cavitationModel (const dictionary &dict, const compressibleTwoPhases &phases) | |
Construct for phases. More... | |
virtual | ~cavitationModel () |
Destructor. More... | |
tmp< volScalarField::Internal > | pSat1 () const |
Return the saturation vapour pressure for phase 1. More... | |
tmp< volScalarField::Internal > | pSat2 () const |
Return the saturation vapour pressure for phase 2. More... | |
Pair< tmp< volScalarField::Internal > > | mDot12Alpha () const |
Return the mass transfer rates of the two phases as coefficients to. More... | |
Pair< tmp< volScalarField::Internal > > | mDot12P () const |
Return the mass transfer rates of the two phases as coefficients to. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from cavitationModel | |
static autoPtr< cavitationModel > | New (const dictionary &dict, const compressibleTwoPhases &phases) |
Protected Member Functions inherited from cavitationModel | |
const volScalarField::Internal & | alphal () const |
Return the liquid density. More... | |
const volScalarField::Internal & | alphav () const |
Return the vapour density. More... | |
const volScalarField::Internal & | rhol () const |
Return the liquid density. More... | |
const volScalarField::Internal & | rhov () const |
Return the vapour density. More... | |
const rhoThermo & | thermol () const |
Return the liquid thermo. More... | |
const rhoThermo & | thermov () const |
Return the vapour thermo. More... | |
tmp< volScalarField::Internal > | pSatl () const |
Return the saturation vapour pressure for the liquid. More... | |
tmp< volScalarField::Internal > | pSatv () const |
Return the saturation vapour pressure for the vapour. More... | |
Protected Attributes inherited from cavitationModel | |
const compressibleTwoPhases & | phases_ |
Phases. More... | |
const bool | liquidIndex_ |
Index of the liquid. More... | |
autoPtr< saturationPressureModel > | saturationModel_ |
The saturation pressure model. More... | |
Saito cavitation model.
Reference:
Saito, Y., Takami, R., Nakamori, I., & Ikohagi, T. (2007). Numerical analysis of unsteady behavior of cloud cavitation around a NACA0015 foil. Computational Mechanics, 40, 85-96.
Usage:
Property | Description | Required | Default value |
---|---|---|---|
liquid | Name of the liquid phase | yes | |
pSat | Saturation vapor pressure | yes | |
Ca | Interfacial area concentration coefficient [1/m] | yes | |
Cv | Vapourisation rate coefficient | yes | |
Cc | Condensation rate coefficient | yes | |
alphaNuc | Nucleation site volume fraction | yes |
Example:
model Saito; liquid liquid; pSat 2300; Ca 0.1; // Interfacial area concentration coefficient [1/m] Cc 1; Cv 1; alphaNuc 0.001;
Saito | ( | const dictionary & | dict, |
const compressibleTwoPhases & | phases | ||
) |
Construct for phases.
Definition at line 52 of file Saito.C.
References Saito::correct().
TypeName | ( | "Saito" | ) |
Runtime type information.
|
virtual |
Return the mass condensation and vaporisation rates as a.
coefficient to multiply alphav for the condensation rate and a coefficient to multiply alphal for the vaporisation rate
Implements cavitationModel.
Definition at line 81 of file Saito.C.
References A, Foam::max(), Foam::min(), and p.
|
virtual |
Return the mass condensation and vaporisation rates as coefficients.
to multiply (p - pSat)
Implements cavitationModel.
Definition at line 111 of file Saito.C.
References A, Foam::max(), Foam::min(), Foam::neg(), p, and Foam::pos0().
|
virtual |
Correct the Saito phaseChange model.
Implements cavitationModel.
Definition at line 140 of file Saito.C.
Referenced by Saito::Saito().
|
virtual |
Read the dictionary and update.
Implements cavitationModel.
Definition at line 144 of file Saito.C.
References dict, and cavitationModel::read().