fvPatch Class Reference

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

Inheritance diagram for fvPatch:

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...
 
 fvPatch (const fvPatch &)
 Disallow default bitwise copy construction. More...
 
virtual ~fvPatch ()
 Destructor. More...
 
const polyPatchpatch () const
 Return the polyPatch. More...
 
virtual const wordname () const
 Return name. More...
 
virtual 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...
 
const objectRegistrydb () const
 Return the local object registry. 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...
 
virtual const scalarFielddeltaCoeffs () const
 Return the face - cell distance coefficient. 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::PatchpatchField (const GeometricField &) const
 Return the corresponding patchField of the named field. More...
 
template<class GeometricField , class Type >
GeometricField::PatchpatchField (GeometricField &) const
 Return the corresponding patchField reference of the named field. More...
 
template<class GeometricField , class Type >
const GeometricField::PatchlookupPatchField (const word &name) const
 Lookup and return the patchField of the named field from the. More...
 
void operator= (const fvPatch &)
 Disallow default bitwise assignment. 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...
 

Friends

class fvBoundaryMesh
 
class surfaceInterpolation
 

Detailed Description

A finiteVolume patch using a polyPatch and a fvBoundaryMesh.

Source files

Definition at line 63 of file fvPatch.H.

Member Typedef Documentation

◆ BoundaryMesh

Definition at line 90 of file fvPatch.H.

Constructor & Destructor Documentation

◆ fvPatch() [1/2]

fvPatch ( const polyPatch p,
const fvBoundaryMesh bm 
)

Construct from polyPatch and fvBoundaryMesh.

Definition at line 46 of file fvPatch.C.

◆ fvPatch() [2/2]

fvPatch ( const fvPatch )

Disallow default bitwise copy construction.

◆ ~fvPatch()

~fvPatch ( )
virtual

Destructor.

Definition at line 55 of file fvPatch.C.

Member Function Documentation

◆ makeWeights()

void makeWeights ( scalarField w) const
protectedvirtual

Make patch weighting factors.

Reimplemented in coupledFvPatch, processorFvPatch, nonConformalProcessorCyclicFvPatch, nonConformalCyclicFvPatch, and cyclicFvPatch.

Definition at line 156 of file fvPatch.C.

◆ initMovePoints()

void initMovePoints ( )
protectedvirtual

Initialise the patches for moving points.

Definition at line 162 of file fvPatch.C.

◆ movePoints()

void movePoints ( )
protectedvirtual

Correct patches after moving points.

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

References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, InfoInFunction, Foam::nl, and fvPatch::patch().

Referenced by fvMeshAdder::add(), fvMesh::addPatch(), and fvMesh::swap().

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

◆ patch()

const polyPatch& patch ( ) const
inline

◆ name()

◆ start()

◆ size()

◆ coupled()

virtual bool coupled ( ) const
inlinevirtual

Return true if this patch is coupled.

Reimplemented in processorFvPatch, nonConformalProcessorCyclicFvPatch, nonConformalCyclicFvPatch, and coupledFvPatch.

Definition at line 163 of file fvPatch.H.

References polyPatch::coupled().

Referenced by Foam::fvc::smooth(), Foam::fvc::spread(), Foam::fvc::sweep(), and fvMeshStitcher::synchronisedBoundaryField().

Here is the call graph for this function:
Here is the caller 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< T >::setSize().

Here is the call graph for this function:

◆ index()

label index ( ) const
inline

Return the index of this patch in the fvBoundaryMesh.

Definition at line 175 of file fvPatch.H.

References patchIdentifier::index().

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::calcPhiYp(), epsilonWallFunctionFvPatchScalarField::calculate(), omegaWallFunctionFvPatchScalarField::calculate(), FvFaceCellWave< Type, TrackingData >::checkCyclic(), fvMeshStitcher::connect(), cyclicFvsPatchField< Type >::cyclicFvsPatchField(), filmFvPatch::deltaCoeffs(), radiationCoupledBase::emissivity(), emptyFvsPatchField< Type >::emptyFvsPatchField(), FvFaceCellWave< Type, TrackingData >::getChangedPatchFaces(), internalFvsPatchField< Type >::internalFvsPatchField(), MRFPatchField::makeAbsolute(), MRFPatchField::makeRelative(), nonConformalProcessorCyclicFvPatch::makeWeights(), FvFaceCellWave< Type, TrackingData >::mergeFaceInfo(), nonConformalCyclicFvsPatchField< Type >::nonConformalCyclicFvsPatchField(), nonConformalErrorFvsPatchField< Type >::nonConformalErrorFvsPatchField(), processorCyclicFvsPatchField< Type >::processorCyclicFvsPatchField(), processorFvsPatchField< Type >::processorFvsPatchField(), symmetryFvsPatchField< Type >::symmetryFvsPatchField(), symmetryPlaneFvsPatchField< Type >::symmetryPlaneFvsPatchField(), fWallFunctionFvPatchScalarField::updateCoeffs(), kLowReWallFunctionFvPatchScalarField::updateCoeffs(), v2WallFunctionFvPatchScalarField::updateCoeffs(), and wedgeFvsPatchField< Type >::wedgeFvsPatchField().

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

◆ boundaryMesh()

◆ db()

const Foam::objectRegistry & db ( ) const

Return the local object registry.

Definition at line 93 of file fvPatch.C.

◆ patchSlice()

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

Slice list to patch.

Definition at line 191 of file fvPatch.H.

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

Referenced by domainDecomposition::procFaceAddressingBf().

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

◆ faceCells()

◆ Cf()

const Foam::vectorField & Cf ( ) const

Return face centres.

Definition at line 105 of file fvPatch.C.

Referenced by nearWallFields::calcAddressing(), and nonConformalProcessorCyclicFvPatch::makeWeights().

Here is the caller graph for this function:

◆ Cn()

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

Return neighbour cell centres.

Definition at line 111 of file fvPatch.C.

References forAll, and tmp< T >::ref().

Here is the call graph for this function:

◆ Sf()

◆ magSf()

const Foam::scalarField & magSf ( ) const

◆ nf()

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

Return face normals.

Definition at line 130 of file fvPatch.C.

Referenced by nearWallFields::calcAddressing(), dynamic::cosTheta(), and nonConformalCoupledFvPatch::makeWeights().

Here is the caller graph for this function:

◆ 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, processorFvPatch, nonConformalProcessorCyclicFvPatch, nonConformalErrorFvPatch, nonConformalCyclicFvPatch, and cyclicFvPatch.

Definition at line 148 of file fvPatch.C.

◆ weights()

const Foam::scalarField & weights ( ) const

Return patch weighting factors.

Definition at line 176 of file fvPatch.C.

◆ deltaCoeffs()

const Foam::scalarField & deltaCoeffs ( ) const
virtual

Return the face - cell distance coefficient.

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

Reimplemented in filmFvPatch.

Definition at line 170 of file fvPatch.C.

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::calcPhiYp(), filmFvPatch::deltaCoeffs(), and filmSurfaceVelocityFvPatchVectorField::updateCoeffs().

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(), forAll, and List< T >::setSize().

Here is the call graph for this function:

◆ patchField() [1/2]

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

Return the corresponding patchField of the named field.

Definition at line 70 of file fvPatchTemplates.C.

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

Here is the call graph for this function:

◆ patchField() [2/2]

GeometricField::Patch & patchField ( GeometricField gf) const

Return the corresponding patchField reference of the named field.

Definition at line 80 of file fvPatchTemplates.C.

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

Here is the call graph for this function:

◆ lookupPatchField()

const GeometricField::Patch & lookupPatchField ( const word name) const

Lookup and return the patchField of the named field from the.

local objectRegistry.

Definition at line 90 of file fvPatchTemplates.C.

References Foam::name().

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::calcPhiYp(), temperatureDependent::cosTheta(), filmSurfaceVelocityFvPatchVectorField::updateCoeffs(), mappedFlowRateVelocityFvPatchVectorField::updateCoeffs(), and coupledTemperatureFvPatchScalarField::updateCoeffs().

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

◆ operator=()

void operator= ( const fvPatch )

Disallow default bitwise assignment.

◆ patchInternalField() [3/3]

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

Definition at line 32 of file fvPatchTemplates.C.

References f(), fvPatch::faceCells(), forAll, tmp< T >::ref(), and fvPatch::size().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ fvBoundaryMesh

friend class fvBoundaryMesh
friend

Definition at line 92 of file fvPatch.H.

◆ surfaceInterpolation

friend class surfaceInterpolation
friend

Definition at line 93 of file fvPatch.H.


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