Saito Class Reference

Saito cavitation model. More...

Inheritance diagram for Saito:
Collaboration diagram for Saito:

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::InternalpSat1 () const
 Return the saturation vapour pressure for phase 1. More...
 
tmp< volScalarField::InternalpSat2 () 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< cavitationModelNew (const dictionary &dict, const compressibleTwoPhases &phases)
 
- Protected Member Functions inherited from cavitationModel
const volScalarField::Internalalphal () const
 Return the liquid density. More...
 
const volScalarField::Internalalphav () const
 Return the vapour density. More...
 
const volScalarField::Internalrhol () const
 Return the liquid density. More...
 
const volScalarField::Internalrhov () const
 Return the vapour density. More...
 
const rhoThermothermol () const
 Return the liquid thermo. More...
 
const rhoThermothermov () const
 Return the vapour thermo. More...
 
tmp< volScalarField::InternalpSatl () const
 Return the saturation vapour pressure for the liquid. More...
 
tmp< volScalarField::InternalpSatv () const
 Return the saturation vapour pressure for the vapour. More...
 
- Protected Attributes inherited from cavitationModel
const compressibleTwoPhasesphases_
 Phases. More...
 
const bool liquidIndex_
 Index of the liquid. More...
 
autoPtr< saturationPressureModelsaturationModel_
 The saturation pressure model. More...
 

Detailed Description

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;
Source files

Definition at line 119 of file Saito.H.

Constructor & Destructor Documentation

◆ Saito()

Saito ( const dictionary dict,
const compressibleTwoPhases phases 
)

Construct for phases.

Definition at line 52 of file Saito.C.

References Saito::correct().

Here is the call graph for this function:

◆ ~Saito()

virtual ~Saito ( )
inlinevirtual

Destructor.

Definition at line 164 of file Saito.H.

Member Function Documentation

◆ TypeName()

TypeName ( "Saito"  )

Runtime type information.

◆ mDotcvAlphal()

Foam::Pair< Foam::tmp< Foam::volScalarField::Internal > > mDotcvAlphal ( ) const
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.

Here is the call graph for this function:

◆ mDotcvP()

Foam::Pair< Foam::tmp< Foam::volScalarField::Internal > > mDotcvP ( ) const
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().

Here is the call graph for this function:

◆ correct()

void correct ( )
virtual

Correct the Saito phaseChange model.

Implements cavitationModel.

Definition at line 140 of file Saito.C.

Referenced by Saito::Saito().

Here is the caller graph for this function:

◆ read()

bool read ( const dictionary dict)
virtual

Read the dictionary and update.

Implements cavitationModel.

Definition at line 144 of file Saito.C.

References dict, and cavitationModel::read().

Here is the call graph for this function:

The documentation for this class was generated from the following files: