Public Member Functions | Protected Member Functions | Private Attributes
immersedBoundaryWallFunctionFvPatchField< Type > Class Template Reference

Boundary condition for passive turbulence variables (U, k, q, R, nut) when using wall functions on an immersed boundary patch. More...

Inheritance diagram for immersedBoundaryWallFunctionFvPatchField< Type >:
Inheritance graph
[legend]
Collaboration diagram for immersedBoundaryWallFunctionFvPatchField< Type >:
Collaboration graph
[legend]

Public Member Functions

 TypeName ("immersedBoundaryWallFunction")
 Runtime type information. More...
 
 immersedBoundaryWallFunctionFvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &)
 Construct from patch and internal field. More...
 
 immersedBoundaryWallFunctionFvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 immersedBoundaryWallFunctionFvPatchField (const immersedBoundaryWallFunctionFvPatchField &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping given. More...
 
 immersedBoundaryWallFunctionFvPatchField (const immersedBoundaryWallFunctionFvPatchField &)
 Construct as copy. More...
 
virtual tmp< fvPatchField< Type > > clone () const
 Construct and return a clone. More...
 
 immersedBoundaryWallFunctionFvPatchField (const immersedBoundaryWallFunctionFvPatchField &, const DimensionedField< Type, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchField< Type > > clone (const DimensionedField< Type, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
virtual ~immersedBoundaryWallFunctionFvPatchField ()
 Destructor. More...
 
Field< Type > & wallValue () const
 Access to value to fix in IB cell. Note non-const access. More...
 
boolListwallMask () const
 Access to indicator on fixed values. Note non-const access. More...
 
- Public Member Functions inherited from immersedBoundaryFvPatchField< Type >
 TypeName (immersedBoundaryFvPatch::typeName_())
 Runtime type information. More...
 
 immersedBoundaryFvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &)
 Construct from patch and internal field. More...
 
 immersedBoundaryFvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 immersedBoundaryFvPatchField (const immersedBoundaryFvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping given immersedBoundaryFvPatchField. More...
 
 immersedBoundaryFvPatchField (const immersedBoundaryFvPatchField< Type > &)
 Construct as copy. More...
 
 immersedBoundaryFvPatchField (const immersedBoundaryFvPatchField< Type > &, const DimensionedField< Type, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual ~immersedBoundaryFvPatchField ()
 Destructor. More...
 
const immersedBoundaryFvPatchibPatch () const
 Return reference to immersed boundary patch. More...
 
virtual Field< Type > & refValue ()
 Return access to reference value. More...
 
virtual const Field< Type > & refValue () const
 Return reference value. More...
 
virtual Field< Type > & refGrad ()
 Return access to reference gradient. More...
 
virtual const Field< Type > & refGrad () const
 Return reference gradient. More...
 
virtual bool fixesValue () const
 Return true if this patch field fixes a value. More...
 
SwitchfixesValue ()
 Return access to fixes value switch. More...
 
const Field< Type > & ibValue () const
 Return fields on intersection points interacting with the IB. More...
 
const Field< Type > & ibGrad () const
 Return IB surface-normal gradient. More...
 
tmp< Field< Type > > ibCellValue () const
 Return IB cell field. More...
 
tmp< Field< Type > > ibSamplingPointValue () const
 Return IB sample point field. More...
 
tmp< Field< Type > > triValue () const
 Return fields on triangular faces. More...
 
tmp< Field< Type > > triGrad () const
 Return IB surface-normal gradient. More...
 
virtual void autoMap (const fvPatchFieldMapper &)
 Map (and resize as needed) from self given a mapping object. More...
 
virtual void rmap (const fvPatchField< Type > &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField. More...
 
void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void initEvaluate (const Pstream::commsTypes commsType=Pstream::blocking)
 Initialise the evaluation of the patch field. More...
 
virtual void evaluate (const Pstream::commsTypes commsType=Pstream::blocking)
 Evaluate the patch field. More...
 
virtual tmp< Field< Type > > valueInternalCoeffs (const tmp< scalarField > &) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< Type > > valueBoundaryCoeffs (const tmp< scalarField > &) const
 Return the matrix source coefficients corresponding to the. More...
 
tmp< Field< Type > > gradientInternalCoeffs () const
 Return the matrix diagonal coefficients corresponding to the. More...
 
tmp< Field< Type > > gradientBoundaryCoeffs () const
 Return the matrix source coefficients corresponding to the. More...
 
virtual void manipulateMatrix (fvMatrix< Type > &matrix)
 Manipulate matrix. More...
 
virtual void write (Ostream &) const
 Write. More...
 
- Public Member Functions inherited from fvPatchField< Type >
 TypeName ("fvPatchField")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< Type > &ptf, const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const dictionary &dict),(p, iF, dict))
 
 fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &)
 Construct from patch and internal field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const word &patchType)
 Construct from patch and internal field and patch type. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &)
 Construct from patch and internal field and patch field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &, const bool valueRequired=false)
 Construct from patch, internal field and dictionary. More...
 
 fvPatchField (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping the given fvPatchField onto a new patch. More...
 
 fvPatchField (const fvPatchField< Type > &)
 Construct as copy. More...
 
 fvPatchField (const fvPatchField< Type > &, const DimensionedField< Type, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchField< Type > > clone () const
 Construct and return a clone. More...
 
Foam::tmp< Foam::fvPatchField< Type > > NewCalculatedType (const fvPatch &p)
 
Foam::tmp< Foam::fvPatchField< Type > > NewCalculatedType (const fvPatchField< Type2 > &pf)
 
virtual ~fvPatchField ()
 Destructor. More...
 
const objectRegistrydb () const
 Return local objectRegistry. More...
 
const fvPatchpatch () const
 Return patch. More...
 
const DimensionedField< Type, volMesh > & dimensionedInternalField () const
 Return dimensioned internal field reference. More...
 
const Field< Type > & internalField () const
 Return internal field reference. More...
 
const wordpatchType () const
 Optional patch type. More...
 
wordpatchType ()
 Optional patch type. More...
 
virtual bool coupled () const
 Return true if this patch field is coupled. More...
 
bool updated () const
 Return true if the boundary condition has already been updated. More...
 
bool manipulatedMatrix () const
 Return true if the matrix has already been manipulated. More...
 
virtual tmp< Field< Type > > snGrad () const
 Return patch-normal gradient. More...
 
virtual tmp< Field< Type > > snGrad (const scalarField &deltaCoeffs) const
 Return patch-normal gradient for coupled-patches. More...
 
virtual void updateCoeffs (const scalarField &weights)
 Update the coefficients associated with the patch field. More...
 
virtual tmp< Field< Type > > patchInternalField () const
 Return internal field next to patch as patch field. More...
 
virtual void patchInternalField (Field< Type > &) const
 Return internal field next to patch as patch field. More...
 
virtual tmp< Field< Type > > patchNeighbourField () const
 Return patchField on the opposite patch of a coupled patch. More...
 
virtual tmp< Field< Type > > gradientInternalCoeffs () const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< Type > > gradientInternalCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< Type > > gradientBoundaryCoeffs () const
 Return the matrix source coefficients corresponding to the. More...
 
virtual tmp< Field< Type > > gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual void manipulateMatrix (fvMatrix< Type > &matrix, const scalarField &weights)
 Manipulate matrix with given weights. More...
 
void writeEntryIfDifferent (Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2) const
 Helper function to write the keyword and entry only if the. More...
 
void check (const fvPatchField< Type > &) const
 Check fvPatchField<Type> against given fvPatchField<Type> More...
 
virtual void operator= (const UList< Type > &)
 
virtual void operator= (const fvPatchField< Type > &)
 
virtual void operator= (const Type &)
 
virtual void operator+= (const fvPatchField< Type > &)
 
virtual void operator+= (const Field< Type > &)
 
virtual void operator+= (const Type &)
 
virtual void operator-= (const fvPatchField< Type > &)
 
virtual void operator-= (const Field< Type > &)
 
virtual void operator-= (const Type &)
 
virtual void operator*= (const fvPatchField< scalar > &)
 
virtual void operator*= (const Field< scalar > &)
 
virtual void operator*= (const scalar)
 
virtual void operator/= (const fvPatchField< scalar > &)
 
virtual void operator/= (const Field< scalar > &)
 
virtual void operator/= (const scalar)
 
virtual void operator== (const fvPatchField< Type > &)
 
virtual void operator== (const Field< Type > &)
 
virtual void operator== (const Type &)
 

Protected Member Functions

virtual void setIbCellValues (const Field< Type > &) const
 Set IB cell values: contains data manipulation. More...
 
- Protected Member Functions inherited from immersedBoundaryFvPatchField< Type >
bool motionUpdateRequired () const
 Does immersed boundary patch field need motion update? More...
 
virtual void motionUpdate () const
 Execute immersed boundary patch field motion update. More...
 

Private Attributes

Field< Type > wallValue_
 Value to fix in IB cell. More...
 
boolList wallMask_
 Indicator on values to fix. More...
 

Additional Inherited Members

- Public Types inherited from fvPatchField< Type >
typedef fvPatch Patch
 
- Static Public Member Functions inherited from fvPatchField< Type >
static tmp< fvPatchField< Type > > New (const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< Type > > New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< Type, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< Type > > New (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
 Return a pointer to a new patchField created on freestore from. More...
 
static tmp< fvPatchField< Type > > New (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
 Return a pointer to a new patchField created on freestore. More...
 
static tmp< fvPatchField< Type > > NewCalculatedType (const fvPatch &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static tmp< fvPatchField< Type > > NewCalculatedType (const fvPatchField< Type2 > &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static const wordcalculatedType ()
 Return the type of the calculated for of fvPatchField. More...
 
- Static Public Attributes inherited from immersedBoundaryFvPatchField< Type >
static const debug::optimisationSwitch nBcIter_
 Number of boundary condition update iterations. More...
 
static const debug::tolerancesSwitch bcTolerance_
 Boundary condition update tolerance. More...
 
- Static Public Attributes inherited from fvPatchField< Type >
static int disallowGenericFvPatchField
 Debug switch to disallow the use of genericFvPatchField. More...
 

Detailed Description

template<class Type>
class Foam::incompressible::RASModels::immersedBoundaryWallFunctionFvPatchField< Type >

Boundary condition for passive turbulence variables (U, k, q, R, nut) when using wall functions on an immersed boundary patch.

Apart form standard immersed boundary capability, the patch field allows only some values to be fixed, as defined by the wall function calculation.

This is handled using the ibWallValue and ibWallMask arrays, corresponding to the number of ibCells

The implementation of wall functions on an immersed boundary will fix the values of k, epsilon, nut and tangential velocity in cells within the log-law layer. For the cells in the laminar sublayer, k and epsilon (omega) will be calculated using a zero gradient condition, nut will be set to zero and the tangential component of the velocity will be corrected.

Author Hrvoje Jasak, Wikki Ltd. All rights reserved

Source files

Definition at line 71 of file immersedBoundaryWallFunctionFvPatchField.H.

Constructor & Destructor Documentation

◆ immersedBoundaryWallFunctionFvPatchField() [1/5]

Construct from patch and internal field.

Definition at line 99 of file immersedBoundaryWallFunctionFvPatchField.C.

Referenced by immersedBoundaryWallFunctionFvPatchField< Type >::clone().

Here is the caller graph for this function:

◆ immersedBoundaryWallFunctionFvPatchField() [2/5]

immersedBoundaryWallFunctionFvPatchField ( const fvPatch p,
const DimensionedField< Type, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 113 of file immersedBoundaryWallFunctionFvPatchField.C.

◆ immersedBoundaryWallFunctionFvPatchField() [3/5]

immersedBoundaryWallFunctionFvPatchField ( const immersedBoundaryWallFunctionFvPatchField< Type > &  ptf,
const fvPatch p,
const DimensionedField< Type, volMesh > &  iF,
const fvPatchFieldMapper mapper 
)

Construct by mapping given.

immersedBoundaryWallFunctionFvPatchField onto a new patch

Definition at line 128 of file immersedBoundaryWallFunctionFvPatchField.C.

◆ immersedBoundaryWallFunctionFvPatchField() [4/5]

Construct as copy.

Definition at line 144 of file immersedBoundaryWallFunctionFvPatchField.C.

◆ immersedBoundaryWallFunctionFvPatchField() [5/5]

Construct as copy setting internal field reference.

Definition at line 157 of file immersedBoundaryWallFunctionFvPatchField.C.

◆ ~immersedBoundaryWallFunctionFvPatchField()

virtual ~immersedBoundaryWallFunctionFvPatchField ( )
inlinevirtual

Destructor.

Definition at line 162 of file immersedBoundaryWallFunctionFvPatchField.H.

Member Function Documentation

◆ setIbCellValues()

void setIbCellValues ( const Field< Type > &  ibcValues) const
protectedvirtual

Set IB cell values: contains data manipulation.

Reimplemented from immersedBoundaryFvPatchField< Type >.

Definition at line 46 of file immersedBoundaryWallFunctionFvPatchField.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, internalField(), Foam::nl, immersedBoundaryFvPatchField< Type >::setIbCellValues(), and List::size().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "immersedBoundaryWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchField<Type> > clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented from immersedBoundaryFvPatchField< Type >.

Definition at line 133 of file immersedBoundaryWallFunctionFvPatchField.H.

References immersedBoundaryWallFunctionFvPatchField< Type >::immersedBoundaryWallFunctionFvPatchField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchField<Type> > clone ( const DimensionedField< Type, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Reimplemented from immersedBoundaryFvPatchField< Type >.

Definition at line 150 of file immersedBoundaryWallFunctionFvPatchField.H.

References immersedBoundaryWallFunctionFvPatchField< Type >::immersedBoundaryWallFunctionFvPatchField().

Here is the call graph for this function:

◆ wallValue()

Foam::Field< Type > & wallValue

Access to value to fix in IB cell. Note non-const access.

Definition at line 171 of file immersedBoundaryWallFunctionFvPatchField.C.

References fvPatch::boundaryMesh(), immersedBoundaryFvPatchField< Type >::ibPatch(), fvBoundaryMesh::mesh(), polyMesh::moving(), and immersedBoundaryFvPatch::movingIb().

Here is the call graph for this function:

◆ wallMask()

Foam::boolList & wallMask

Access to indicator on fixed values. Note non-const access.

Definition at line 200 of file immersedBoundaryWallFunctionFvPatchField.C.

References fvPatch::boundaryMesh(), immersedBoundaryFvPatchField< Type >::ibPatch(), fvBoundaryMesh::mesh(), polyMesh::moving(), and immersedBoundaryFvPatch::movingIb().

Here is the call graph for this function:

Field Documentation

◆ wallValue_

Field<Type> wallValue_
mutableprivate

Value to fix in IB cell.

Definition at line 78 of file immersedBoundaryWallFunctionFvPatchField.H.

◆ wallMask_

boolList wallMask_
mutableprivate

Indicator on values to fix.

Definition at line 81 of file immersedBoundaryWallFunctionFvPatchField.H.


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