Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Friends
fvPatch Class Reference

A finiteVolume patch using a polyPatch and a fvBoundaryMesh. More...

Inheritance diagram for fvPatch:
Inheritance graph
[legend]
Collaboration diagram for fvPatch:
Collaboration graph
[legend]

Public Types

typedef fvBoundaryMesh BoundaryMesh
 

Public Member Functions

 TypeName (polyPatch::typeName_())
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
 
 fvPatch (const polyPatch &, const fvBoundaryMesh &)
 Construct from polyPatch and fvBoundaryMesh. More...
 
virtual ~fvPatch ()
 Destructor. More...
 
const polyPatchpatch () const
 Return the polyPatch. More...
 
const wordname () const
 Return name. More...
 
label start () const
 Return start label of this patch in the polyMesh face list. More...
 
virtual label size () const
 Return size. More...
 
virtual bool coupled () const
 Return true if this patch is coupled. More...
 
label index () const
 Return the index of this patch in the fvBoundaryMesh. More...
 
const fvBoundaryMeshboundaryMesh () const
 Return boundaryMesh reference. More...
 
template<class T >
const List< T >::subList patchSlice (const List< T > &l) const
 Slice list to patch. More...
 
virtual const labelUListfaceCells () const
 Return faceCells. More...
 
const vectorFieldCf () const
 Return face centres. More...
 
tmp< vectorFieldCn () const
 Return neighbour cell centres. More...
 
const vectorFieldSf () const
 Return face area vectors. More...
 
const scalarFieldmagSf () const
 Return face area magnitudes. More...
 
tmp< vectorFieldnf () const
 Return face normals. More...
 
virtual tmp< vectorFielddelta () const
 Return cell-centre to face-centre vector. More...
 
const scalarFieldweights () const
 Return patch weighting factors. More...
 
const scalarFielddeltaCoeffs () const
 Return the face - cell distance coeffient. More...
 
template<class Type >
tmp< Field< Type > > patchInternalField (const UList< Type > &) const
 Return given internal field next to patch as patch field. More...
 
template<class Type >
void patchInternalField (const UList< Type > &, Field< Type > &) const
 Return given internal field next to patch as patch field. More...
 
template<class GeometricField , class Type >
const GeometricField::PatchFieldTypepatchField (const GeometricField &) const
 Return the corresponding patchField of the named field. More...
 
template<class GeometricField , class Type >
const GeometricField::PatchFieldTypelookupPatchField (const word &name, const GeometricField *=NULL, const Type *=NULL) const
 Lookup and return the patchField of the named field from the. More...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > patchInternalField (const UList< Type > &f) const
 

Static Public Member Functions

static autoPtr< fvPatchNew (const polyPatch &, const fvBoundaryMesh &)
 Return a pointer to a new patch created on freestore from polyPatch. More...
 
static bool constraintType (const word &pt)
 Return true if the given type is a constraint type. More...
 
static wordList constraintTypes ()
 Return a list of all the constraint patch types. More...
 

Protected Member Functions

virtual void makeWeights (scalarField &) const
 Make patch weighting factors. More...
 
virtual void initMovePoints ()
 Initialise the patches for moving points. More...
 
virtual void movePoints ()
 Correct patches after moving points. More...
 

Private Member Functions

 fvPatch (const fvPatch &)
 Disallow construct as copy. More...
 
void operator= (const fvPatch &)
 Disallow assignment. More...
 

Private Attributes

const polyPatchpolyPatch_
 Reference to the underlying polyPatch. More...
 
const fvBoundaryMeshboundaryMesh_
 Reference to boundary mesh. More...
 

Friends

class fvBoundaryMesh
 
class surfaceInterpolation
 

Detailed Description

A finiteVolume patch using a polyPatch and a fvBoundaryMesh.

Source files

Definition at line 61 of file fvPatch.H.

Member Typedef Documentation

◆ BoundaryMesh

Definition at line 97 of file fvPatch.H.

Constructor & Destructor Documentation

◆ fvPatch() [1/2]

fvPatch ( const fvPatch )
private

Disallow construct as copy.

◆ fvPatch() [2/2]

fvPatch ( const polyPatch p,
const fvBoundaryMesh bm 
)

Construct from polyPatch and fvBoundaryMesh.

Definition at line 46 of file fvPatch.C.

◆ ~fvPatch()

~fvPatch ( )
virtual

Destructor.

Definition at line 55 of file fvPatch.C.

Member Function Documentation

◆ operator=()

void operator= ( const fvPatch )
private

Disallow assignment.

◆ makeWeights()

void makeWeights ( scalarField w) const
protectedvirtual

Make patch weighting factors.

Reimplemented in coupledFvPatch, cyclicACMIFvPatch, cyclicFvPatch, cyclicAMIFvPatch, and processorFvPatch.

Definition at line 150 of file fvPatch.C.

References w().

Referenced by cyclicAMIFvPatch::makeWeights(), and cyclicACMIFvPatch::makeWeights().

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

◆ initMovePoints()

void initMovePoints ( )
protectedvirtual

Initialise the patches for moving points.

Reimplemented in immersedBoundaryFvPatch.

Definition at line 156 of file fvPatch.C.

◆ movePoints()

void movePoints ( )
protectedvirtual

Correct patches after moving points.

Reimplemented in immersedBoundaryFvPatch.

Definition at line 160 of file fvPatch.C.

◆ TypeName()

TypeName ( polyPatch::typeName_()  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
fvPatch  ,
polyPatch  ,
(const polyPatch &patch, const fvBoundaryMesh &bm)  ,
(patch, bm)   
)

◆ New()

Foam::autoPtr< Foam::fvPatch > New ( const polyPatch patch,
const fvBoundaryMesh bm 
)
static

Return a pointer to a new patch created on freestore from polyPatch.

Definition at line 33 of file fvPatchNew.C.

References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::Info, and Foam::nl.

Referenced by fvBoundaryMesh::addPatches(), and meshRefinement::appendPatch().

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

◆ patch()

const polyPatch& patch ( ) const
inline

◆ name()

const word& name ( ) const
inline

◆ start()

label start ( ) const
inline

◆ size()

virtual label size ( ) const
inlinevirtual

◆ coupled()

virtual bool coupled ( ) const
inlinevirtual

Return true if this patch is coupled.

Reimplemented in cyclicACMIFvPatch, cyclicAMIFvPatch, regionCoupledWallFvPatch, regionCoupledFvPatch, processorFvPatch, and coupledFvPatch.

Definition at line 167 of file fvPatch.H.

References polyPatch::coupled(), and fvPatch::polyPatch_.

Here is the call graph for this function:

◆ constraintType()

bool constraintType ( const word pt)
static

Return true if the given type is a constraint type.

Definition at line 61 of file fvPatch.C.

◆ constraintTypes()

Foam::wordList constraintTypes ( )
static

Return a list of all the constraint patch types.

Definition at line 67 of file fvPatch.C.

References List::setSize().

Here is the call graph for this function:

◆ index()

label index ( ) const
inline

◆ boundaryMesh()

const fvBoundaryMesh& boundaryMesh ( ) const
inline

◆ patchSlice()

const List<T>::subList patchSlice ( const List< T > &  l) const
inline

Slice list to patch.

Definition at line 192 of file fvPatch.H.

References fvPatch::size(), and fvPatch::start().

Referenced by fvPatchMapper::calcAddressing().

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

◆ faceCells()

const Foam::labelUList & faceCells ( ) const
virtual

Return faceCells.

Reimplemented in regionCoupledWallFvPatch, regionCoupledFvPatch, coupledFvPatch, and emptyFvPatch.

Definition at line 93 of file fvPatch.C.

Referenced by curvatureSeparation::calcCosAngle(), nearWallDist::calculate(), epsilonLowReWallFunctionFvPatchScalarField::calculate(), epsilonWallFunctionFvPatchScalarField::calculate(), omegaWallFunctionFvPatchScalarField::calculate(), epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields(), omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields(), inversePointDistanceDiffusivity::correct(), inverseFaceDistanceDiffusivity::correct(), contactAngleForce::correct(), SSG< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), pointLinear< Type >::correction(), coupledFvPatch::faceCells(), regionCoupledFvPatch::faceCells(), regionCoupledWallFvPatch::faceCells(), epsilonWallFunctionFvPatchScalarField::manipulateMatrix(), omegaWallFunctionFvPatchScalarField::manipulateMatrix(), energyRegionCoupledFvPatchScalarField::patchNeighbourField(), energyRegionCoupledFvPatchScalarField::patchNeighbourTemperatureField(), energyJumpAMIFvPatchScalarField::updateCoeffs(), energyJumpFvPatchScalarField::updateCoeffs(), fWallFunctionFvPatchScalarField::updateCoeffs(), v2WallFunctionFvPatchScalarField::updateCoeffs(), kLowReWallFunctionFvPatchScalarField::updateCoeffs(), epsilonWallFunctionFvPatchScalarField::updateCoeffs(), and omegaWallFunctionFvPatchScalarField::updateCoeffs().

Here is the caller graph for this function:

◆ Cf()

const Foam::vectorField & Cf ( ) const

Return face centres.

Definition at line 99 of file fvPatch.C.

References boundaryMesh::mesh().

Referenced by meshToMesh0::interpolate(), and main().

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

◆ Cn()

Foam::tmp< Foam::vectorField > Cn ( ) const

Return neighbour cell centres.

Definition at line 105 of file fvPatch.C.

References forAll, and boundaryMesh::mesh().

Here is the call graph for this function:

◆ Sf()

const Foam::vectorField & Sf ( ) const

Return face area vectors.

Definition at line 130 of file fvPatch.C.

References boundaryMesh::mesh().

Referenced by cyclicACMIFvPatch::updateAreas(), activeBaffleVelocityFvPatchVectorField::updateCoeffs(), and activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs().

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

◆ magSf()

const Foam::scalarField & magSf ( ) const

◆ nf()

Foam::tmp< Foam::vectorField > nf ( ) const

◆ delta()

Foam::tmp< Foam::vectorField > delta ( ) const
virtual

Return cell-centre to face-centre vector.

except for coupled patches for which the cell-centre to coupled-cell-centre vector is returned

Reimplemented in coupledFvPatch, cyclicACMIFvPatch, cyclicAMIFvPatch, processorFvPatch, and cyclicFvPatch.

Definition at line 142 of file fvPatch.C.

Referenced by cyclicACMIFvPatch::delta(), cyclicACMIFvPatch::makeWeights(), and energyRegionCoupledFvPatchScalarField::weights().

Here is the caller graph for this function:

◆ weights()

const Foam::scalarField & weights ( ) const

Return patch weighting factors.

Definition at line 170 of file fvPatch.C.

References boundaryMesh::mesh().

Here is the call graph for this function:

◆ deltaCoeffs()

const Foam::scalarField & deltaCoeffs ( ) const

Return the face - cell distance coeffient.

except for coupled patches for which the cell-centre to coupled-cell-centre distance coeffient is returned

Definition at line 164 of file fvPatch.C.

References boundaryMesh::mesh().

Referenced by contactAngleForce::correct(), totalFlowRateAdvectiveDiffusiveFvPatchScalarField::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs(), turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs(), and humidityTemperatureCoupledMixedFvPatchScalarField::updateCoeffs().

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

◆ patchInternalField() [1/3]

tmp<Field<Type> > patchInternalField ( const UList< Type > &  ) const

Return given internal field next to patch as patch field.

◆ patchInternalField() [2/3]

void patchInternalField ( const UList< Type > &  f,
Field< Type > &  pif 
) const

Return given internal field next to patch as patch field.

Definition at line 52 of file fvPatchTemplates.C.

References f(), and forAll.

Here is the call graph for this function:

◆ patchField()

const GeometricField::PatchFieldType & patchField ( const GeometricField gf) const

Return the corresponding patchField of the named field.

Definition at line 70 of file fvPatchTemplates.C.

References GeometricField::boundaryField().

Referenced by porousBafflePressureFvPatchField< Type >::updateCoeffs().

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

◆ lookupPatchField()

const GeometricField::PatchFieldType & lookupPatchField ( const word name,
const GeometricField = NULL,
const Type *  = NULL 
) const

◆ patchInternalField() [3/3]

Foam::tmp<Foam::Field<Type> > patchInternalField ( const UList< Type > &  f) const

Definition at line 32 of file fvPatchTemplates.C.

References f(), and forAll.

Here is the call graph for this function:

Friends And Related Function Documentation

◆ fvBoundaryMesh

friend class fvBoundaryMesh
friend

Definition at line 99 of file fvPatch.H.

◆ surfaceInterpolation

friend class surfaceInterpolation
friend

Definition at line 100 of file fvPatch.H.

Field Documentation

◆ polyPatch_

const polyPatch& polyPatch_
private

Reference to the underlying polyPatch.

Definition at line 66 of file fvPatch.H.

Referenced by fvPatch::coupled(), fvPatch::index(), fvPatch::name(), fvPatch::patch(), fvPatch::size(), and fvPatch::start().

◆ boundaryMesh_

const fvBoundaryMesh& boundaryMesh_
private

Reference to boundary mesh.

Definition at line 69 of file fvPatch.H.

Referenced by fvPatch::boundaryMesh().


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