alphatWallBoilingWallFunctionFvPatchScalarField Class Reference

A thermal wall function for simulation of subcooled nucleate wall boiling with runtime selectable submodels for: More...

Inheritance diagram for alphatWallBoilingWallFunctionFvPatchScalarField:
Collaboration diagram for alphatWallBoilingWallFunctionFvPatchScalarField:

Classes

struct  boilingLiquidProperties
 
struct  properties
 

Public Types

enum  phaseType { vapourPhase , vaporPhase , liquidPhase }
 Enumeration listing the possible operational modes. More...
 

Public Member Functions

 TypeName ("compressible::alphatWallBoilingWallFunction")
 Runtime type information. More...
 
 alphatWallBoilingWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 alphatWallBoilingWallFunctionFvPatchScalarField (const alphatWallBoilingWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fieldMapper &)
 Construct by mapping given. More...
 
 alphatWallBoilingWallFunctionFvPatchScalarField (const alphatWallBoilingWallFunctionFvPatchScalarField &)=delete
 Disallow copy without setting internal field reference. More...
 
 alphatWallBoilingWallFunctionFvPatchScalarField (const alphatWallBoilingWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Copy constructor setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
const scalarFieldwetFraction () const
 Return the wall liquid fraction field [-]. More...
 
const scalarFielddDeparture () const
 Return the departure diameter field [m]. More...
 
const scalarFieldfDeparture () const
 Return the departure frequency field [Hz]. More...
 
const scalarFieldnucleationSiteDensity () const
 Return the nucleation site density field [1/m^2]. More...
 
const scalarFieldqQuenching () const
 Return the quenching surface heat flux field [W/m^2]. More...
 
const scalarFieldqEvaporative () const
 Return the evaporative surface heat flux field [W/m^2]. More...
 
virtual const scalarFielddmdtf () const
 Return the rate of phase change. More...
 
virtual void map (const fvPatchScalarField &, const fieldMapper &)
 Map the given fvPatchField onto this fvPatchField. More...
 
virtual void reset (const fvPatchScalarField &)
 Reset the fvPatchField to the given fvPatchField. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void write (Ostream &) const
 Write. More...
 
- Public Member Functions inherited from alphatPhaseChangeWallFunctionBase
 TypeName ("compressible::alphatPhaseChangeWallFunctionBase")
 Runtime type information. More...
 
 alphatPhaseChangeWallFunctionBase ()
 Construct null. More...
 
 alphatPhaseChangeWallFunctionBase (const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const dictionary &)
 Construct from a patch, an internal field and a dictionary. More...
 
virtual ~alphatPhaseChangeWallFunctionBase ()
 
bool activeInterface (const phaseInterface &) const
 Is there phase change mass transfer for this interface? More...
 

Static Public Attributes

static const NamedEnum< phaseType, 3 > phaseTypeNames_
 Heat source type names. More...
 

Additional Inherited Members

- Protected Attributes inherited from alphatPhaseChangeWallFunctionBase
const word phaseName_
 Name of the phase. More...
 
const word otherPhaseName_
 Name of the other phase. More...
 

Detailed Description

A thermal wall function for simulation of subcooled nucleate wall boiling with runtime selectable submodels for:

- wall heat flux partitioning model

  • nucleation site density
  • bubble departure frequency
  • bubble departure diameter

Implements a version of the well-known RPI wall boiling model (Kurul & Podowski, 1991). The model implementation based on implementation described in Peltola et al. (2019) and is similar to the model described by Peltola & Pättikangas (2012). The present implementation includes simplified support for presence of non-volatile components in addition to a single volatile component.

References:

    Kurul, N., & Podowski, M.Z. (1991).
    On the modeling of multidimensional effects in boiling channels.
    ANS. Proc. National Heat Transfer Con. Minneapolis, Minnesota, USA,
    1991.
    ISBN: 0-89448-162-1, pp. 30-40.
    Peltola, J., Pättikangas, T., Bainbridge, W., Lehnigk, R., Schlegel, F.
    (2019).
    On Development and validation of subcooled nucleate boiling models for
    OpenFOAM Foundation Release.
    NURETH-18 Conference Proceedings, Portland, Oregon, United States, 2019.
    Peltola, J., & Pättikangas, T.J.H. (2012).
    Development and validation of a boiling model for OpenFOAM multiphase
    solver.
    CFD4NRS-4 Conference Proceedings, Daejeon, Korea, 2012.
    paper 59.
Usage
Property Description Required Default value
phaseType 'vapour' or 'liquid' yes
useLiquidTemperatureWallFunction Use wall function to calculate liquid temperature? no yes
tolerance solution tolerance no rootSmall
Prt turbulent Prandtl number no 0.85

if phaseType 'vapour':

partitioningModel yes

if phaseType 'liquid':

partitioningModel yes
nucleationSiteModel yes
departureDiameterModel yes
departureFrequenctModel yes
bubbleWaitingTimeRatio no 0.8

NOTE: Runtime selectable submodels may require model specific entries

Example usage:

    hotWall
    {
        type            compressible::alphatWallBoilingWallFunction;
        phaseType       liquid;
        dmdt            uniform 0;
        Prt             0.85;
        partitioningModel
        {
            type        Lavieville;
            alphaCrit   0.2;
        }
        nucleationSiteModel
        {
            type        LemmertChawla;
        }
        departureDiameterModel
        {
            type        TolubinskiKostanchuk;
        }
        departureFrequencyModel
        {
            type        Cole;
        }
        value           uniform 0.01;
Source files

Definition at line 200 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

Member Enumeration Documentation

◆ phaseType

enum phaseType

Enumeration listing the possible operational modes.

Enumerator
vapourPhase 
vaporPhase 
liquidPhase 

Definition at line 210 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ alphatWallBoilingWallFunctionFvPatchScalarField() [1/4]

◆ alphatWallBoilingWallFunctionFvPatchScalarField() [2/4]

Construct by mapping given.

alphatWallBoilingWallFunctionFvPatchScalarField onto a new patch

Definition at line 708 of file alphatWallBoilingWallFunctionFvPatchScalarField.C.

◆ alphatWallBoilingWallFunctionFvPatchScalarField() [3/4]

Disallow copy without setting internal field reference.

◆ alphatWallBoilingWallFunctionFvPatchScalarField() [4/4]

Copy constructor setting internal field reference.

Definition at line 742 of file alphatWallBoilingWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ TypeName()

TypeName ( "compressible::alphatWallBoilingWallFunction"  )

Runtime type information.

◆ clone()

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Definition at line 390 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

References alphatWallBoilingWallFunctionFvPatchScalarField::alphatWallBoilingWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ wetFraction()

const scalarField& wetFraction ( ) const
inline

Return the wall liquid fraction field [-].

Definition at line 405 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ dDeparture()

const scalarField& dDeparture ( ) const
inline

Return the departure diameter field [m].

Definition at line 411 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ fDeparture()

const scalarField& fDeparture ( ) const
inline

Return the departure frequency field [Hz].

Definition at line 417 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ nucleationSiteDensity()

const scalarField& nucleationSiteDensity ( ) const
inline

Return the nucleation site density field [1/m^2].

Definition at line 423 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ qQuenching()

const scalarField& qQuenching ( ) const
inline

Return the quenching surface heat flux field [W/m^2].

Definition at line 429 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ qEvaporative()

const scalarField& qEvaporative ( ) const
inline

Return the evaporative surface heat flux field [W/m^2].

Definition at line 435 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ dmdtf()

virtual const scalarField& dmdtf ( ) const
inlinevirtual

Return the rate of phase change.

Implements alphatPhaseChangeWallFunctionBase.

Definition at line 441 of file alphatWallBoilingWallFunctionFvPatchScalarField.H.

◆ map()

void map ( const fvPatchScalarField ptf,
const fieldMapper mapper 
)
virtual

Map the given fvPatchField onto this fvPatchField.

Definition at line 776 of file alphatWallBoilingWallFunctionFvPatchScalarField.C.

◆ reset()

void reset ( const fvPatchScalarField ptf)
virtual

Reset the fvPatchField to the given fvPatchField.

Used for mesh to mesh mapping

Definition at line 797 of file alphatWallBoilingWallFunctionFvPatchScalarField.C.

References Field< Type >::reset().

Here is the call graph for this function:

◆ updateCoeffs()

◆ write()

Member Data Documentation

◆ phaseTypeNames_

const NamedEnum<phaseType, 3> phaseTypeNames_
static

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