curvatureSeparation Class Reference

Curvature film separation model. More...

Inheritance diagram for curvatureSeparation:
[legend]
Collaboration diagram for curvatureSeparation:
[legend]

Public Member Functions

 TypeName ("curvatureSeparation")
 Runtime type information. More...
 
 curvatureSeparation (liquidFilmBase &film, const dictionary &dict)
 Construct from surface film model. More...
 
virtual ~curvatureSeparation ()=default
 Destructor. More...
 
virtual void correct (scalarField &availableMass, scalarField &massToInject, scalarField &diameterToInject)
 Correct. More...
 
- Public Member Functions inherited from injectionModel
 TypeName ("injectionModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, injectionModel, dictionary,(liquidFilmBase &film, const dictionary &dict),(film, dict))
 
 injectionModel (liquidFilmBase &film)
 Construct null. More...
 
 injectionModel (const word &modelType, liquidFilmBase &film, const dictionary &dict)
 Construct from type name, dictionary and surface film model. More...
 
virtual ~injectionModel ()
 Destructor. More...
 
virtual void correct (scalarField &availableMass, scalarField &massToInject, scalarField &diameterToInject)=0
 Correct. More...
 
virtual scalar injectedMassTotal () const
 Return the total mass injected. More...
 
virtual void patchInjectedMassTotals (scalar &patchMasses) const
 
- Public Member Functions inherited from filmSubModelBase
 filmSubModelBase (liquidFilmBase &film)
 Construct null. More...
 
 filmSubModelBase (liquidFilmBase &film, const dictionary &dict, const word &baseName, const word &modelType, const word &dictExt="Coeffs")
 Construct from film film without name. More...
 
 filmSubModelBase (const word &modelName, liquidFilmBase &film, const dictionary &dict, const word &baseName, const word &modelType)
 Construct from film film with name. More...
 
virtual ~filmSubModelBase ()
 Destructor. More...
 
virtual bool writeTime () const
 Flag to indicate when to write a property. More...
 
const liquidFilmBasefilm () const
 Return const access to the film surface film model. More...
 
liquidFilmBasefilm ()
 Return the reference to the film surface film model. More...
 
template<class FilmType >
const FilmType & filmType () const
 
- Public Member Functions inherited from subModelBase
 subModelBase (dictionary &properties)
 Construct null. More...
 
 subModelBase (dictionary &properties, const dictionary &dict, const word &baseName, const word &modelType, const word &dictExt="Coeffs")
 Construct from components without name. More...
 
 subModelBase (const word &modelName, dictionary &properties, const dictionary &dict, const word &baseName, const word &modelType)
 Construct from components with name. More...
 
 subModelBase (const subModelBase &smb)
 Construct as copy. More...
 
virtual ~subModelBase ()
 Destructor. More...
 
const wordmodelName () const
 Return const access to the name of the sub-model. More...
 
const dictionarydict () const
 Return const access to the cloud dictionary. More...
 
const wordbaseName () const
 Return const access to the base name of the sub-model. More...
 
const wordmodelType () const
 Return const access to the sub-model type. More...
 
const dictionarycoeffDict () const
 Return const access to the coefficients dictionary. More...
 
const dictionaryproperties () const
 Return const access to the properties dictionary. More...
 
virtual bool defaultCoeffs (const bool printMsg) const
 Returns true if defaultCoeffs is true and outputs on printMsg. More...
 
virtual bool active () const
 Return the model 'active' status - default active = true. More...
 
virtual void cacheFields (const bool store)
 Cache dependent sub-model fields. More...
 
virtual bool writeTime () const
 Flag to indicate when to write a property. More...
 
virtual fileName localPath () const
 Output directory. More...
 
template<class Type >
Type getBaseProperty (const word &entryName, const Type &defaultValue=Type(Zero)) const
 Retrieve generic property from the base model. More...
 
template<class Type >
void getBaseProperty (const word &entryName, Type &value) const
 Retrieve generic property from the base model. More...
 
template<class Type >
void setBaseProperty (const word &entryName, const Type &value)
 Add generic property to the base model. More...
 
bool getModelDict (const word &entryName, dictionary &dict) const
 Retrieve dictionary, return true if set. More...
 
template<class Type >
bool getModelProperty (const word &entryName, Type &value) const
 Retrieve generic property from the sub-model. More...
 
template<class Type >
Type getModelProperty (const word &entryName, const Type &defaultValue=Type(Zero)) const
 Retrieve generic property from the sub-model. More...
 
template<class Type >
void setModelProperty (const word &entryName, const Type &value)
 Add generic property to the sub-model. More...
 
virtual void write (Ostream &os) const
 Write. More...
 

Protected Member Functions

tmp< areaScalarFieldcalcInvR1 (const areaVectorField &U, const scalarField &calcCosAngle) const
 Calculate local (inverse) radius of curvature. More...
 
tmp< scalarFieldcalcCosAngle (const edgeScalarField &phi) const
 
- Protected Member Functions inherited from injectionModel
void addToInjectedMass (const scalar dMass)
 Add to injected mass. More...
 
void correct ()
 Correct. More...
 
- Protected Member Functions inherited from subModelBase
bool inLine () const
 Flag to indicate whether data is/was read in-line. More...
 

Protected Attributes

areaTensorField gradNHat_
 Gradient of surface normals. More...
 
scalar deltaByR1Min_
 Minimum gravity driven film thickness (non-dimensionalised delta/R1) More...
 
scalar definedPatchRadii_
 
scalar magG_
 Magnitude of gravity vector. More...
 
vector gHat_
 Direction of gravity vector. More...
 
scalar fThreshold_
 Threshold force for separation. More...
 
scalar minInvR1_
 Minimum inv R1 for separation. More...
 
- Protected Attributes inherited from filmSubModelBase
liquidFilmBasefilmModel_
 Reference to the film surface film model. More...
 
- Protected Attributes inherited from subModelBase
const word modelName_
 Name of the sub-model. More...
 
dictionaryproperties_
 Reference to properties dictionary e.g. for restart. More...
 
const dictionary dict_
 Copy of dictionary used during construction. More...
 
const word baseName_
 Name of the sub-model base class. More...
 
const word modelType_
 Type of the sub-model. More...
 
const dictionary coeffDict_
 Coefficients dictionary. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from injectionModel
static autoPtr< injectionModelNew (liquidFilmBase &film, const dictionary &dict, const word &mdoelType)
 Return a reference to the selected injection model. More...
 

Detailed Description

Curvature film separation model.

Assesses film curvature via the mesh geometry and calculates a force balance of the form:

F_sum = F_inertial + F_body + F_surface

If F_sum < 0, the film separates. Similarly, if F_sum > 0 the film will remain attached.

Based on description given by Owen and D. J. Ryley. The flow of thin liquid films around corners. International Journal of Multiphase Flow, 11(1):51-62, 1985.

Source files

Definition at line 68 of file curvatureSeparation.H.

Constructor & Destructor Documentation

◆ curvatureSeparation()

curvatureSeparation ( liquidFilmBase film,
const dictionary dict 
)

Construct from surface film model.

Definition at line 157 of file curvatureSeparation.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, filmSubModelBase::film(), liquidFilmBase::g(), curvatureSeparation::gHat_, curvatureSeparation::magG_, and dimensioned< Type >::value().

Here is the call graph for this function:

◆ ~curvatureSeparation()

virtual ~curvatureSeparation ( )
virtualdefault

Destructor.

Member Function Documentation

◆ calcInvR1()

tmp< areaScalarField > calcInvR1 ( const areaVectorField U,
const scalarField calcCosAngle 
) const
protected

Calculate local (inverse) radius of curvature.

Definition at line 56 of file curvatureSeparation.C.

References curvatureSeparation::definedPatchRadii_, Foam::dimVelocity, e, forAll, curvatureSeparation::gradNHat_, Foam::mag(), Foam::max(), tmp< T >::ref(), and U.

Referenced by curvatureSeparation::calcCosAngle(), and curvatureSeparation::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calcCosAngle()

tmp< scalarField > calcCosAngle ( const edgeScalarField phi) const
protected

Calculate the cosine of the angle between gravity vector and cell out flow direction

Definition at line 94 of file curvatureSeparation.C.

References curvatureSeparation::calcInvR1(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimless, Foam::dimVelocity, filmSubModelBase::film(), forAll, curvatureSeparation::gHat_, Foam::mag(), Foam::max(), mesh, Foam::min(), primitiveMesh::nFaces(), IOobject::NO_READ, phi, Foam::pos(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), regionFaModel::regionMesh(), UList< T >::size(), fvMesh::time(), timeName, U, liquidFilmBase::Uf(), regIOobject::write(), TimeState::writeTime(), and Foam::Zero.

Referenced by curvatureSeparation::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "curvatureSeparation"  )

Runtime type information.

◆ correct()

Member Data Documentation

◆ gradNHat_

areaTensorField gradNHat_
protected

Gradient of surface normals.

Definition at line 86 of file curvatureSeparation.H.

Referenced by curvatureSeparation::calcInvR1().

◆ deltaByR1Min_

scalar deltaByR1Min_
protected

Minimum gravity driven film thickness (non-dimensionalised delta/R1)

Definition at line 89 of file curvatureSeparation.H.

Referenced by curvatureSeparation::correct().

◆ definedPatchRadii_

scalar definedPatchRadii_
protected

List of radii for patches - if patch not defined, radius calculated based on mesh geometry

Definition at line 93 of file curvatureSeparation.H.

Referenced by curvatureSeparation::calcInvR1().

◆ magG_

scalar magG_
protected

Magnitude of gravity vector.

Definition at line 96 of file curvatureSeparation.H.

Referenced by curvatureSeparation::correct(), and curvatureSeparation::curvatureSeparation().

◆ gHat_

vector gHat_
protected

Direction of gravity vector.

Definition at line 99 of file curvatureSeparation.H.

Referenced by curvatureSeparation::calcCosAngle(), and curvatureSeparation::curvatureSeparation().

◆ fThreshold_

scalar fThreshold_
protected

Threshold force for separation.

Definition at line 102 of file curvatureSeparation.H.

Referenced by curvatureSeparation::correct().

◆ minInvR1_

scalar minInvR1_
protected

Minimum inv R1 for separation.

Definition at line 105 of file curvatureSeparation.H.

Referenced by curvatureSeparation::correct().


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