A finiteVolume patch using a polyPatch and a fvBoundaryMesh. More...
#include <fvPatch.H>

Public Types | |
| typedef fvBoundaryMesh | BoundaryMesh |
| The boundary type associated with the patch. | |
Public Member Functions | |
| virtual void | makeWeights (scalarField &) const |
| Make patch weighting factors. | |
| virtual void | makeDeltaCoeffs (scalarField &) const |
| Correct patch deltaCoeffs. | |
| virtual void | makeNonOrthoDeltaCoeffs (scalarField &) const |
| Correct patch non-ortho deltaCoeffs. | |
| virtual void | makeNonOrthoCorrVectors (vectorField &) const |
| Correct patch non-ortho correction vectors. | |
| virtual void | initMovePoints () |
| Initialise the patches for moving points. | |
| virtual void | movePoints () |
| Correct patches after moving points. | |
| TypeName (polyPatch::typeName_()) | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm)) | |
| fvPatch (const polyPatch &, const fvBoundaryMesh &) | |
| Construct from polyPatch and fvBoundaryMesh. | |
| virtual | ~fvPatch () |
| Destructor. | |
| const polyPatch & | patch () const noexcept |
| Return the polyPatch. | |
| virtual const word & | name () const |
| Return name. | |
| label | index () const noexcept |
| The index of this patch in the boundary mesh. | |
| label | start () const noexcept |
| The patch start within the polyMesh face list. | |
| label | offset () const noexcept |
| The offset of this patch within the boundary face list. | |
| virtual label | size () const |
| Patch size is the number of faces, but can be overloaded. | |
| virtual bool | coupled () const |
| Return true if this patch is coupled. | |
| const fvBoundaryMesh & | boundaryMesh () const noexcept |
| Return boundaryMesh reference. | |
| template<class T> | |
| const List< T >::subList | patchSlice (const UList< T > &values) const |
| This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size. | |
| template<class T> | |
| const List< T >::subList | boundarySlice (const UList< T > &values) const |
| This patch slice from the list of boundary values, which has size mesh::nBoundaryFaces(), using the virtual patch size. | |
| virtual const labelUList & | faceCells () const |
| Return faceCells. | |
| const vectorField & | Cf () const |
| Return face centres. | |
| tmp< vectorField > | Cn () const |
| Return neighbour cell centres. | |
| const vectorField & | Sf () const |
| Return face area vectors, like the fvMesh::Sf() method. | |
| const scalarField & | magSf () const |
| Return face area magnitudes, like the fvMesh::magSf() method. | |
| tmp< vectorField > | unitSf () const |
| Return face unit normals, like the fvMesh::unitSf() method. Same as nf(). | |
| tmp< vectorField > | nf () const |
| Return face unit normals, like the fvMesh::unitSf() method Same as unitSf(). | |
| virtual tmp< vectorField > | delta () const |
| Return cell-centre to face-centre vector except for coupled patches for which the cell-centre to coupled-cell-centre vector is returned. | |
| const scalarField & | weights () const |
| Return patch weighting factors. | |
| const 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. | |
| template<class Type> | |
| void | patchInternalField (const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const |
| Extract internal field next to patch using specified addressing. | |
| template<class Type> | |
| void | patchInternalField (const UList< Type > &internalData, UList< Type > &pfld) const |
| Extract internal field next to patch as patch field using faceCells() mapping. | |
| 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. | |
| template<class GeometricField, class AnyType = bool> | |
| const GeometricField::Patch & | patchField (const GeometricField &gf) const |
| Return the patch field of the GeometricField corresponding to this patch. | |
| template<class GeometricField, class AnyType = bool> | |
| const GeometricField::Patch & | lookupPatchField (const word &name) const |
| Lookup the named field from the local registry and return the patch field corresponding to this patch. | |
| template<class GeometricField> | |
| const GeometricField::Patch * | cfindPatchField (const word &name) const |
| Find the named field (if any) from the local registry and return the patch field corresponding to this patch. | |
| template<class Type> | |
| Foam::tmp< Foam::Field< Type > > | patchInternalField (const UList< Type > &internalData) const |
Static Public Member Functions | |
| static autoPtr< fvPatch > | New (const polyPatch &, const fvBoundaryMesh &) |
| Return a pointer to a new patch created on freestore from polyPatch. | |
| static const fvPatch & | lookupPatch (const polyPatch &p) |
| Lookup the polyPatch index on corresponding fvMesh. | |
| static bool | constraintType (const word &patchType) |
| Return true if the given type is a constraint type. | |
| static wordList | constraintTypes () |
| Return a list of all the constraint patch types. | |
Friends | |
| class | fvBoundaryMesh |
| class | surfaceInterpolation |
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
| typedef fvBoundaryMesh BoundaryMesh |
| fvPatch | ( | const polyPatch & | p, |
| const fvBoundaryMesh & | bm ) |
Construct from polyPatch and fvBoundaryMesh.
Definition at line 59 of file fvPatch.C.
References fvBoundaryMesh, and p.

|
virtual |
Make patch weighting factors.
Reimplemented in coupledFvPatch, cyclicACMIFvPatch, cyclicAMIFvPatch, cyclicFvPatch, and processorFvPatch.
Definition at line 157 of file fvPatch.C.
Referenced by cyclicACMIFvPatch::makeWeights(), and cyclicAMIFvPatch::makeWeights().

|
virtual |
Correct patch deltaCoeffs.
Reimplemented in cyclicAMIFvPatch, and mappedVariableThicknessWallFvPatch.
|
virtual |
Correct patch non-ortho deltaCoeffs.
Reimplemented in cyclicAMIFvPatch.
|
virtual |
Correct patch non-ortho correction vectors.
Reimplemented in cyclicAMIFvPatch.
|
virtual |
|
virtual |
Correct patches after moving points.
Reimplemented in cyclicACMIFvPatch, and cyclicAMIFvPatch.
| TypeName | ( | polyPatch::typeName_() | ) |
Runtime type information.
|
static |
Return a pointer to a new patch created on freestore from polyPatch.
Definition at line 27 of file fvPatchNew.C.
References DebugInFunction, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInLookup, fvBoundaryMesh, and patch().
Referenced by fvMeshAdder::add(), fvMeshTools::addPatch(), and meshRefinement::appendPatch().


|
static |
Lookup the polyPatch index on corresponding fvMesh.
Definition at line 42 of file fvPatch.C.
References fvMesh::boundary(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::isA(), Foam::nl, and p.

|
inlinenoexcept |
Return the polyPatch.
Definition at line 202 of file fvPatch.H.
References Foam::noexcept.
Referenced by reconstructedDistanceFunction::constructRDF(), coupledFvPatch::coupledFvPatch(), cyclicACMIFvPatch::cyclicACMIFvPatch(), cyclicAMIFvPatch::cyclicAMIFvPatch(), cyclicFvPatch::cyclicFvPatch(), cyclicSlipFvPatch::cyclicSlipFvPatch(), declareRunTimeSelectionTable(), emptyFvPatch::emptyFvPatch(), genericFvPatch::genericFvPatch(), fvMeshSubset::interpolate(), mappedVariableThicknessWallFvPatch::makeDeltaCoeffs(), mappedFvPatch::mappedFvPatch(), mappedVariableThicknessWallFvPatch::mappedVariableThicknessWallFvPatch(), mappedWallFvPatch::mappedWallFvPatch(), inverseDistance::markBoundaries(), trackingInverseDistance::markBoundaries(), cyclicACMIFvPatch::movePoints(), cyclicAMIFvPatch::movePoints(), New(), nonuniformTransformCyclicFvPatch::nonuniformTransformCyclicFvPatch(), oversetFvPatch::oversetFvPatch(), processorCyclicFvPatch::processorCyclicFvPatch(), processorFvPatch::processorFvPatch(), cyclicACMIFvPatch::resetPatchAreas(), symmetryFvPatch::symmetryFvPatch(), symmetryPlaneFvPatch::symmetryPlaneFvPatch(), cyclicACMIFvPatch::TypeName(), emptyFvPatch::TypeName(), wedgeFvPatch::TypeName(), activeBaffleVelocityFvPatchVectorField::updateCoeffs(), activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs(), outletMappedUniformInletHeatAdditionFvPatchField::updateCoeffs(), wallFvPatch::wallFvPatch(), wedgeFvPatch::wedgeFvPatch(), and regionSizeDistribution::write().

|
inlinevirtual |
Return name.
Definition at line 210 of file fvPatch.H.
Referenced by advectiveFvPatchField< Type >::advectiveFvPatchField(), fieldExtents::calcFieldExtents(), genericFvPatchField< Type >::genericFvPatchField(), genericFvPatchField< Type >::genericFvPatchField(), genericFvsPatchField< Type >::genericFvsPatchField(), genericFvsPatchField< Type >::genericFvsPatchField(), meshToMesh0::interpolate(), oversetFvPatchField< Type >::manipulateMatrix(), patchInternalField(), fvMeshDistribute::testField(), uniformMixedFvPatchField< Type >::uniformMixedFvPatchField(), cyclicACMIFvPatch::updateAreas(), thermalBaffle1DFvPatchScalarField< solidType >::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), humidityTemperatureCoupledMixedFvPatchScalarField::updateCoeffs(), and mappedFlowRateFvPatchVectorField::updateCoeffs().

|
inlinenoexcept |
The index of this patch in the boundary mesh.
Definition at line 218 of file fvPatch.H.
References Foam::noexcept.
Referenced by Cf(), cfindPatchField(), cyclicACMIFvsPatchField< Type >::cyclicACMIFvsPatchField(), cyclicACMIFvsPatchField< Type >::cyclicACMIFvsPatchField(), cyclicAMIFvsPatchField< Type >::cyclicAMIFvsPatchField(), cyclicAMIFvsPatchField< Type >::cyclicAMIFvsPatchField(), cyclicFvsPatchField< Type >::cyclicFvsPatchField(), cyclicFvsPatchField< Type >::cyclicFvsPatchField(), deltaCoeffs(), emptyFvsPatchField< Type >::emptyFvsPatchField(), emptyFvsPatchField< Type >::emptyFvsPatchField(), lookupPatchField(), magSf(), cyclicAMIFvPatch::movePoints(), patchField(), processorCyclicFvsPatchField< Type >::processorCyclicFvsPatchField(), processorCyclicFvsPatchField< Type >::processorCyclicFvsPatchField(), processorFvsPatchField< Type >::processorFvsPatchField(), processorFvsPatchField< Type >::processorFvsPatchField(), Sf(), symmetryFvsPatchField< Type >::symmetryFvsPatchField(), symmetryFvsPatchField< Type >::symmetryFvsPatchField(), symmetryPlaneFvsPatchField< Type >::symmetryPlaneFvsPatchField(), symmetryPlaneFvsPatchField< Type >::symmetryPlaneFvsPatchField(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), wedgeFvsPatchField< Type >::wedgeFvsPatchField(), wedgeFvsPatchField< Type >::wedgeFvsPatchField(), and weights().

|
inlinenoexcept |
The patch start within the polyMesh face list.
Definition at line 226 of file fvPatch.H.
References Foam::noexcept.
Referenced by CentredFitSnGradData< Polynomial >::calcFit(), LeastSquaresGrad< Type, Stencil >::calcGrad(), extendedCellToFaceStencil::collectData(), extendedFaceToCellStencil::collectData(), adjointBoundaryCondition< Type >::computePatchGrad(), inverseFaceDistanceDiffusivity::correct(), inversePointDistanceDiffusivity::correct(), pointLinear< Type >::correction(), emptyFvPatch::emptyFvPatch(), cellDecomposer::interpolate(), fvMeshSubset::interpolate(), fvMeshSubset::interpolate(), dynamicRefineFvMesh::mapFields(), fvMeshAdder::MapSurfaceField(), Foam::oversetAdjustPhi(), patchSlice(), oversetFvPatchField< Type >::storeFringeCoefficients(), fvMeshDistribute::testField(), mappedFixedInternalValueFvPatchField< Type >::updateCoeffs(), mappedVelocityFluxFixedValueFvPatchField::updateCoeffs(), extendedCellToFaceStencil::weightedSum(), and extendedUpwindCellToFaceStencil::weightedSum().

|
inlinenoexcept |
The offset of this patch within the boundary face list.
Definition at line 234 of file fvPatch.H.
References Foam::noexcept.
Referenced by boundarySlice().

|
inlinevirtual |
Patch size is the number of faces, but can be overloaded.
Reimplemented in emptyFvPatch.
Definition at line 242 of file fvPatch.H.
Referenced by boundarySlice(), cyclicACMIFvPatch::coupled(), cyclicAMIFvPatch::delta(), cyclicACMIFvPatch::interfaceInternalField(), cyclicAMIFvPatch::interfaceInternalField(), cyclicFvPatch::interfaceInternalField(), oversetFvPatch::interfaceInternalField(), processorFvPatch::interfaceInternalField(), processorFvPatch::internalFieldTransfer(), cellDecomposer::interpolate(), fvMeshSubset::interpolate(), fvMeshSubset::interpolate(), cyclicAMIFvPatch::makeWeights(), patchInternalField(), patchInternalField(), and patchSlice().

|
inlinevirtual |
Return true if this patch is coupled.
Reimplemented in coupledFvPatch, cyclicACMIFvPatch, cyclicAMIFvPatch, and processorFvPatch.
|
static |
Return true if the given type is a constraint type.
Definition at line 74 of file fvPatch.C.
References UList< T >::found().
Referenced by constraintTypes(), cellVolumeWeight::findHoles(), inverseDistance::findHoles(), inverseDistance::markBoundaries(), trackingInverseDistance::markBoundaries(), and cellVolumeWeight::markPatchCells().


|
static |
Return a list of all the constraint patch types.
Definition at line 85 of file fvPatch.C.
References constraintType(), forAllConstIters, and List< T >::resize().

|
inlinenoexcept |
Return boundaryMesh reference.
Definition at line 268 of file fvPatch.H.
References fvBoundaryMesh, and Foam::noexcept.
Referenced by Cf(), cfindPatchField(), Cn(), cyclicAMIFvPatch::coupled(), deltaCoeffs(), emptyFvPatch::emptyFvPatch(), lookupPatchField(), magSf(), waWallFunctionFvPatchScalarField::manipulateMatrix(), cyclicACMIFvPatch::neighbFvPatch(), cyclicAMIFvPatch::neighbFvPatch(), cyclicFvPatch::neighbFvPatch(), cyclicACMIFvPatch::neighbPatch(), cyclicAMIFvPatch::neighbPatch(), cyclicFvPatch::neighbPatch(), cyclicACMIFvPatch::nonOverlapPatch(), Sf(), outletMappedUniformInletHeatAdditionFvPatchField::updateCoeffs(), and weights().


|
virtual |
Return faceCells.
Reimplemented in coupledFvPatch, emptyFvPatch, and oversetFvPatch.
Definition at line 107 of file fvPatch.C.
Referenced by epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields(), omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields(), inverseFaceDistanceDiffusivity::correct(), inversePointDistanceDiffusivity::correct(), LRR< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), contactAngleForce::correct(), ThermalPhaseChangePhaseSystem< BasePhaseSystem >::correctInterfaceThermo(), pointLinear< Type >::correction(), coupledFvPatch::faceCells(), cellVolumeWeight::findHoles(), inverseDistance::findHoles(), waWallFunctionFvPatchScalarField::manipulateMatrix(), inverseDistance::markBoundaries(), trackingInverseDistance::markBoundaries(), cellVolumeWeight::markPatchCells(), Foam::oversetAdjustPhi(), patchInternalField(), patchInternalField(), and oversetFvPatchField< Type >::storeFringeCoefficients().

| const Foam::vectorField & Cf | ( | ) | const |
Return face centres.
Definition at line 113 of file fvPatch.C.
References boundaryMesh(), and index().
Referenced by fieldExtents::calcFieldExtents(), coupledFvPatch::delta(), delta(), meshToMesh0::interpolate(), cyclicAMIFvPatch::movePoints(), and cyclicACMIFvPatch::resetPatchAreas().


| Foam::tmp< Foam::vectorField > Cn | ( | ) | const |
Return neighbour cell centres.
Definition at line 119 of file fvPatch.C.
References boundaryMesh(), mesh, and patchInternalField().
Referenced by coupledFvPatch::delta(), and delta().


| const Foam::vectorField & Sf | ( | ) | const |
Return face area vectors, like the fvMesh::Sf() method.
Definition at line 125 of file fvPatch.C.
References boundaryMesh(), and index().
Referenced by cyclicAMIFvPatch::movePoints(), nf(), cyclicACMIFvPatch::resetPatchAreas(), unitSf(), activeBaffleVelocityFvPatchVectorField::updateCoeffs(), and activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs().


| const Foam::scalarField & magSf | ( | ) | const |
Return face area magnitudes, like the fvMesh::magSf() method.
Definition at line 131 of file fvPatch.C.
References boundaryMesh(), and index().
Referenced by cyclicACMIFvPatch::movePoints(), cyclicAMIFvPatch::movePoints(), nf(), cyclicACMIFvPatch::resetPatchAreas(), unitSf(), activeBaffleVelocityFvPatchVectorField::updateCoeffs(), activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs(), outletMappedUniformInletFvPatchField< Type >::updateCoeffs(), and outletMappedUniformInletHeatAdditionFvPatchField::updateCoeffs().


| Foam::tmp< Foam::vectorField > unitSf | ( | ) | const |
| Foam::tmp< Foam::vectorField > nf | ( | ) | const |
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition at line 143 of file fvPatch.C.
Referenced by ReynoldsAnalogy::Cf(), delta(), cyclicACMIFvPatch::makeWeights(), cyclicAMIFvPatch::makeWeights(), cyclicFvPatch::makeWeights(), and processorFvPatch::makeWeights().


|
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, cyclicFvPatch, and processorFvPatch.
Definition at line 149 of file fvPatch.C.
References Cf(), Cn(), and nf().

| const Foam::scalarField & weights | ( | ) | const |
Return patch weighting factors.
Definition at line 189 of file fvPatch.C.
References boundaryMesh(), and index().

| 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(), and index().
Referenced by contactAngleForce::correct(), turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), and humidityTemperatureCoupledMixedFvPatchScalarField::updateCoeffs().


|
inline |
Extract internal field next to patch using specified addressing.
| internalData | The internal field to extract from | |
| addressing | Addressing from patch into internal field | |
| [out] | pfld | The extracted patch field. Should normally be sized according to the patch size(), which can be smaller than the addressing size |
Definition at line 27 of file fvPatchTemplates.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, FOAM_UNLIKELY, Foam::nl, size(), and UList< T >::size().
Referenced by Cn(), processorFvPatch::initInternalFieldTransfer(), cyclicACMIFvPatch::interfaceInternalField(), cyclicACMIFvPatch::interfaceInternalField(), cyclicAMIFvPatch::interfaceInternalField(), cyclicAMIFvPatch::interfaceInternalField(), cyclicFvPatch::interfaceInternalField(), cyclicFvPatch::interfaceInternalField(), oversetFvPatch::interfaceInternalField(), oversetFvPatch::interfaceInternalField(), processorFvPatch::interfaceInternalField(), processorFvPatch::interfaceInternalField(), patchInternalField(), and patchInternalField().


| void patchInternalField | ( | const UList< Type > & | internalData, |
| UList< Type > & | pfld ) const |
Extract internal field next to patch as patch field using faceCells() mapping.
| internalData | The internal field to extract from | |
| [out] | pfld | The extracted patch field. Should normally be sized according to the patch size(), which can be smaller than the faceCells() size |
Definition at line 61 of file fvPatchTemplates.C.
References faceCells(), and patchInternalField().

|
nodiscard |
Return given internal field next to patch as patch field using faceCells() mapping.
| internalData | The internal field to extract from |
References name().

| const GeometricField::Patch & patchField | ( | const GeometricField & | gf | ) | const |
Return the patch field of the GeometricField corresponding to this patch.
Definition at line 89 of file fvPatchTemplates.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), and index().

| const GeometricField::Patch & lookupPatchField | ( | const word & | name | ) | const |
Lookup the named field from the local registry and return the patch field corresponding to this patch.
Definition at line 28 of file fvPatchFvMeshTemplates.C.
References boundaryMesh(), index(), fvBoundaryMesh::mesh(), and fvMesh::thisDb().
Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::phiY(), turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs(), filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs(), humidityTemperatureCoupledMixedFvPatchScalarField::updateCoeffs(), and mappedFlowRateFvPatchVectorField::updateCoeffs().


| const GeometricField::Patch * cfindPatchField | ( | const word & | name | ) | const |
Find the named field (if any) from the local registry and return the patch field corresponding to this patch.
Definition at line 39 of file fvPatchFvMeshTemplates.C.
References boundaryMesh(), and index().

| Foam::tmp< Foam::Field< Type > > patchInternalField | ( | const UList< Type > & | internalData | ) | const |
Definition at line 77 of file fvPatchTemplates.C.
References faceCells(), Foam::New(), patchInternalField(), and size().

|
friend |
Definition at line 140 of file fvPatch.H.
References fvBoundaryMesh.
Referenced by boundaryMesh(), coupledFvPatch::coupledFvPatch(), declareRunTimeSelectionTable(), emptyFvPatch::emptyFvPatch(), fvBoundaryMesh, fvPatch(), genericFvPatch::genericFvPatch(), mappedFvPatch::mappedFvPatch(), New(), oversetFvPatch::oversetFvPatch(), symmetryFvPatch::symmetryFvPatch(), symmetryPlaneFvPatch::symmetryPlaneFvPatch(), emptyFvPatch::TypeName(), wedgeFvPatch::TypeName(), wallFvPatch::wallFvPatch(), and wedgeFvPatch::wedgeFvPatch().
|
friend |
Definition at line 141 of file fvPatch.H.
References surfaceInterpolation.
Referenced by surfaceInterpolation.