Public Types | Public Member Functions | Static Public Member Functions | Friends | List of all members
fvPatch Class Reference

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

Inheritance diagram for fvPatch:
Inheritance graph
[legend]

Public Types

typedef fvBoundaryMesh BoundaryMesh
 The boundary type associated with the patch. More...
 

Public Member Functions

virtual void makeWeights (scalarField &) const
 Make patch weighting factors. More...
 
virtual void makeDeltaCoeffs (scalarField &) const
 Correct patch deltaCoeffs. More...
 
virtual void makeNonOrthoDeltaCoeffs (scalarField &) const
 Correct patch non-ortho deltaCoeffs. More...
 
virtual void makeNonOrthoCorrVectors (vectorField &) const
 Correct patch non-ortho correction vectors. More...
 
virtual void initMovePoints ()
 Initialise the patches for moving points. More...
 
virtual void movePoints ()
 Correct patches after moving points. More...
 
 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 noexcept
 Return the polyPatch. More...
 
virtual const wordname () const
 Return name. More...
 
label index () const noexcept
 The index of this patch in the boundary mesh. More...
 
label start () const noexcept
 The patch start within the polyMesh face list. More...
 
virtual label size () const
 Patch size is the number of faces, but can be overloaded. More...
 
virtual bool coupled () const
 Return true if this patch is coupled. More...
 
const fvBoundaryMeshboundaryMesh () const noexcept
 Return boundaryMesh reference. More...
 
template<class T >
const List< T >::subList patchSlice (const List< T > &values) const
 This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size. More...
 
template<class T >
const List< T >::subList boundarySlice (const List< T > &values) const
 This patch slice from the list of boundary values, which has size mesh::nBoundaryFaces(), using the virtual patch size. 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, like the fvMesh::Sf() method. More...
 
const scalarFieldmagSf () const
 Return face area magnitudes, like the fvMesh::magSf() method. More...
 
tmp< vectorFieldunitSf () const
 Return face unit normals, like the fvMesh::unitSf() method. Same as nf(). More...
 
tmp< vectorFieldnf () const
 Return face unit normals, like the fvMesh::unitSf() method Same as unitSf(). More...
 
virtual tmp< vectorFielddelta () const
 Return cell-centre to face-centre vector except for coupled patches for which the cell-centre to coupled-cell-centre vector is returned. More...
 
const scalarFieldweights () const
 Return patch weighting factors. More...
 
const scalarFielddeltaCoeffs () const
 Return the face - cell distance coefficient except for coupled patches for which the cell-centre to coupled-cell-centre distance coefficient is returned. More...
 
template<class Type >
void patchInternalField (const UList< Type > &internalData, const labelUList &addressing, Field< Type > &pfld) const
 Extract internal field next to patch using specified addressing. More...
 
template<class Type >
void patchInternalField (const UList< Type > &internalData, Field< Type > &pfld) const
 Extract internal field next to patch as patch field using faceCells() mapping. More...
 
template<class Type >
tmp< Field< Type > > patchInternalField (const UList< Type > &internalData) const
 Return given internal field next to patch as patch field using faceCells() mapping. More...
 
template<class GeometricField , class AnyType = bool>
const GeometricField::PatchpatchField (const GeometricField &gf) const
 Return the patch field of the GeometricField corresponding to this patch. More...
 
template<class GeometricField , class AnyType = bool>
const GeometricField::PatchlookupPatchField (const word &name, const GeometricField *=nullptr, const AnyType *=nullptr) const
 Lookup the named field from the local registry and return the patch field corresponding to this patch. More...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > patchInternalField (const UList< Type > &internalData) 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 const fvPatchlookupPatch (const polyPatch &p)
 Lookup the polyPatch index on corresponding fvMesh. More...
 
static bool constraintType (const word &patchType)
 Return true if the given type is a constraint type. More...
 
static wordList constraintTypes ()
 Return a list of all the constraint patch types. More...
 

Friends

class fvBoundaryMesh
 
class surfaceInterpolation
 

Detailed Description

A finiteVolume patch using a polyPatch and a fvBoundaryMesh.

Source files

Definition at line 70 of file fvPatch.H.

Member Typedef Documentation

◆ BoundaryMesh

The boundary type associated with the patch.

Definition at line 138 of file fvPatch.H.

Constructor & Destructor Documentation

◆ fvPatch()

fvPatch ( const polyPatch p,
const fvBoundaryMesh bm 
)

Construct from polyPatch and fvBoundaryMesh.

Definition at line 59 of file fvPatch.C.

◆ ~fvPatch()

~fvPatch ( )
virtual

Destructor.

Definition at line 68 of file fvPatch.C.

Member Function Documentation

◆ makeWeights()

void makeWeights ( scalarField w) const
virtual

Make patch weighting factors.

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

Definition at line 157 of file fvPatch.C.

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

Here is the caller graph for this function:

◆ makeDeltaCoeffs()

void makeDeltaCoeffs ( scalarField w) const
virtual

Correct patch deltaCoeffs.

Reimplemented in cyclicAMIFvPatch, and mappedVariableThicknessWallFvPatch.

Definition at line 163 of file fvPatch.C.

◆ makeNonOrthoDeltaCoeffs()

void makeNonOrthoDeltaCoeffs ( scalarField w) const
virtual

Correct patch non-ortho deltaCoeffs.

Reimplemented in cyclicAMIFvPatch.

Definition at line 167 of file fvPatch.C.

◆ makeNonOrthoCorrVectors()

void makeNonOrthoCorrVectors ( vectorField w) const
virtual

Correct patch non-ortho correction vectors.

Reimplemented in cyclicAMIFvPatch.

Definition at line 171 of file fvPatch.C.

◆ initMovePoints()

void initMovePoints ( )
virtual

Initialise the patches for moving points.

Definition at line 175 of file fvPatch.C.

◆ movePoints()

void movePoints ( )
virtual

Correct patches after moving points.

Reimplemented in cyclicACMIFvPatch, and cyclicAMIFvPatch.

Definition at line 179 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 28 of file fvPatchNew.C.

References DebugInFunction, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInLookup, and Foam::foamVersion::patch.

Referenced by fvMeshAdder::add(), and meshRefinement::appendPatch().

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

◆ lookupPatch()

const Foam::fvPatch & lookupPatch ( const polyPatch p)
static

Lookup the polyPatch index on corresponding fvMesh.

Note
Fatal if the polyPatch is not associated with a fvMesh

Definition at line 42 of file fvPatch.C.

References fvMesh::boundary(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::nl, and p.

Here is the call graph for this function:

◆ patch()

const polyPatch& patch ( ) const
inlinenoexcept

◆ name()

virtual const word& name ( ) const
inlinevirtual

◆ index()

label index ( ) const
inlinenoexcept

The index of this patch in the boundary mesh.

Definition at line 218 of file fvPatch.H.

References patchIdentifier::index().

Referenced by boundaryAdjointContributionIncompressible::laminarDiffusivity(), boundaryAdjointContributionIncompressible::momentumDiffusion(), cyclicAMIFvPatch::movePoints(), boundaryAdjointContributionIncompressible::pab(), boundaryAdjointContributionIncompressible::pb(), boundaryAdjointContributionIncompressible::phiab(), boundaryAdjointContributionIncompressible::phib(), boundaryAdjointContributionIncompressible::pressureSource(), boundaryAdjointContributionIncompressible::tangentVelocitySource(), boundaryAdjointContributionIncompressible::TMVariable1(), boundaryAdjointContributionIncompressible::TMVariable1Diffusion(), boundaryAdjointContributionIncompressible::TMVariable2(), boundaryAdjointContributionIncompressible::TMVariable2Diffusion(), boundaryAdjointContributionIncompressible::turbulentDiffusivity(), boundaryAdjointContributionIncompressible::Uab(), boundaryAdjointContributionIncompressible::Ub(), energyJumpAMIFvPatchScalarField::updateCoeffs(), energyJumpFvPatchScalarField::updateCoeffs(), vibrationShellFvPatchScalarField::updateCoeffs(), totalFlowRateAdvectiveDiffusiveFvPatchScalarField::updateCoeffs(), thermalShellFvPatchScalarField::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), kLowReWallFunctionFvPatchScalarField::updateCoeffs(), velocityFilmShellFvPatchVectorField::updateCoeffs(), boundaryAdjointContributionIncompressible::velocitySource(), and boundaryAdjointContributionIncompressible::wallDistance().

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

◆ start()

label start ( ) const
inlinenoexcept

◆ size()

virtual label size ( ) const
inlinevirtual

◆ coupled()

virtual bool coupled ( ) const
inlinevirtual

Return true if this patch is coupled.

Reimplemented in cyclicAMIFvPatch, cyclicACMIFvPatch, processorFvPatch, and coupledFvPatch.

Definition at line 242 of file fvPatch.H.

References polyPatch::coupled().

Here is the call graph for this function:

◆ constraintType()

bool constraintType ( const word patchType)
static

Return true if the given type is a constraint type.

Definition at line 74 of file fvPatch.C.

References UList< T >::found().

Here is the call graph for this function:

◆ constraintTypes()

Foam::wordList constraintTypes ( )
static

Return a list of all the constraint patch types.

Definition at line 85 of file fvPatch.C.

References forAllConstIters().

Here is the call graph for this function:

◆ boundaryMesh()

const fvBoundaryMesh& boundaryMesh ( ) const
inlinenoexcept

◆ patchSlice()

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

This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size.

Definition at line 271 of file fvPatch.H.

References fvPatch::size(), fvPatch::start(), and Foam::HashTableOps::values().

Referenced by fvFieldDecomposer::reset().

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

◆ boundarySlice()

const List<T>::subList boundarySlice ( const List< T > &  values) const
inline

This patch slice from the list of boundary values, which has size mesh::nBoundaryFaces(), using the virtual patch size.

Definition at line 282 of file fvPatch.H.

References polyPatch::offset(), fvPatch::size(), and Foam::HashTableOps::values().

Here is the call graph for this function:

◆ faceCells()

const Foam::labelUList & faceCells ( ) const
virtual

◆ Cf()

const Foam::vectorField & Cf ( ) const

Return face centres.

Definition at line 113 of file fvPatch.C.

References boundaryMesh::mesh().

Referenced by fieldExtents::calcFieldExtents(), coupledFvPatch::delta(), meshToMesh0::interpolate(), cyclicAMIFvPatch::movePoints(), and cyclicACMIFvPatch::resetPatchAreas().

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 119 of file fvPatch.C.

References mesh.

Referenced by coupledFvPatch::delta().

Here is the caller graph for this function:

◆ Sf()

const Foam::vectorField & Sf ( ) const

Return face area vectors, like the fvMesh::Sf() method.

Definition at line 125 of file fvPatch.C.

References boundaryMesh::mesh().

Referenced by cyclicAMIFvPatch::movePoints(), cyclicACMIFvPatch::resetPatchAreas(), 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

◆ unitSf()

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

Return face unit normals, like the fvMesh::unitSf() method. Same as nf().

Definition at line 137 of file fvPatch.C.

◆ 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 cyclicAMIFvPatch, cyclicACMIFvPatch, processorFvPatch, cyclicFvPatch, and coupledFvPatch.

Definition at line 149 of file fvPatch.C.

◆ weights()

const Foam::scalarField & weights ( ) const

Return patch weighting factors.

Definition at line 189 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 coefficient except for coupled patches for which the cell-centre to coupled-cell-centre distance coefficient is returned.

Definition at line 183 of file fvPatch.C.

References boundaryMesh::mesh().

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

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

◆ patchInternalField() [1/4]

void patchInternalField ( const UList< Type > &  internalData,
const labelUList addressing,
Field< Type > &  pfld 
) const
inline

Extract internal field next to patch using specified addressing.

Parameters
internalDataThe internal field to extract from
addressingAddressing from patch into internal field
[out]pfldThe extracted patch field. It is always resized according to the patch size(), which can be smaller than the addressing size

Definition at line 28 of file fvPatchTemplates.C.

References List< Type >::resize_nocopy().

Here is the call graph for this function:

◆ patchInternalField() [2/4]

void patchInternalField ( const UList< Type > &  internalData,
Field< Type > &  pfld 
) const

Extract internal field next to patch as patch field using faceCells() mapping.

Parameters
internalDataThe internal field to extract from
[out]pfldThe extracted patch field. It is always resized according to the patch size(), which can be smaller than the faceCells() size

Definition at line 47 of file fvPatchTemplates.C.

◆ patchInternalField() [3/4]

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

Return given internal field next to patch as patch field using faceCells() mapping.

Parameters
internalDataThe internal field to extract from

◆ patchField()

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

Return the patch field of the GeometricField corresponding to this patch.

Definition at line 70 of file fvPatchTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField().

Here is the call graph for this function:

◆ lookupPatchField()

const GeometricField::Patch & lookupPatchField ( const word name,
const GeometricField = nullptr,
const AnyType *  = nullptr 
) const

Lookup the named field from the local registry and return the patch field corresponding to this patch.

N.B. The dummy pointer arguments are used if this function is instantiated within a templated function to avoid a bug in gcc.

Definition at line 28 of file fvPatchFvMeshTemplates.C.

References boundaryMesh::mesh(), and Foam::name().

Referenced by enthalpySorptionFvPatchScalarField::patchSource(), semiPermeableBaffleMassFractionFvPatchScalarField::phiY(), totalFlowRateAdvectiveDiffusiveFvPatchScalarField::updateCoeffs(), mappedFlowRateFvPatchVectorField::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), porousBafflePressureFvPatchField::updateCoeffs(), turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs(), and humidityTemperatureCoupledMixedFvPatchScalarField::updateCoeffs().

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

◆ patchInternalField() [4/4]

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

Definition at line 58 of file fvPatchTemplates.C.

References Foam::New().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ fvBoundaryMesh

friend class fvBoundaryMesh
friend

Definition at line 140 of file fvPatch.H.

◆ surfaceInterpolation

friend class surfaceInterpolation
friend

Definition at line 141 of file fvPatch.H.


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