fvMeshSubset Class Reference

Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the subMesh) with mapping lists for points, faces, and cells. More...

Inheritance diagram for fvMeshSubset:
Collaboration diagram for fvMeshSubset:

Public Member Functions

 fvMeshSubset (const fvMesh &baseMesh)
 Construct using the entire mesh (no subset) More...
 
 fvMeshSubset (const fvMesh &baseMesh, const Foam::zero)
 Construct a zero-sized subset mesh, non-processor patches only. More...
 
 fvMeshSubset (const fvMesh &baseMesh, const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Construct for a cell-subset of the given mesh. More...
 
 fvMeshSubset (const fvMesh &baseMesh, const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
 Construct for a cell-subset of the given mesh. More...
 
 fvMeshSubset (const fvMesh &baseMesh, const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Construct for a cell-subset of the given mesh. More...
 
 fvMeshSubset (const fvMesh &baseMesh, const label regioni, const labelUList &regions, const label patchID=-1, const bool syncPar=true)
 Construct for a cell-subset of the given mesh. More...
 
const fvMeshbaseMesh () const noexcept
 Original mesh. More...
 
const fvMeshmesh () const noexcept
 Return baseMesh or subMesh, depending on the current state. More...
 
bool hasSubMesh () const noexcept
 Have subMesh? More...
 
const fvMeshsubMesh () const
 Return reference to subset mesh. More...
 
fvMeshsubMesh ()
 Return reference to subset mesh. More...
 
const labelListpointMap () const
 Return point map. More...
 
const labelListfaceMap () const
 Return face map. More...
 
const labelListfaceFlipMap () const
 Return face map with sign to encode flipped faces. More...
 
const labelListcellMap () const
 Return cell map. More...
 
const labelListpatchMap () const
 Return patch map. More...
 
void clear ()
 Reset subMesh and all maps. More...
 
void reset ()
 Reset subMesh and all maps. Same as clear() More...
 
void reset (const Foam::zero)
 Reset to a zero-sized subset mesh, non-processor patches only. More...
 
void reset (autoPtr< fvMesh > &&subMeshPtr, labelList &&pointMap, labelList &&faceMap, labelList &&cellMap, labelList &&patchMap)
 Reset from components. More...
 
void reset (const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. More...
 
void reset (const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. More...
 
void reset (const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. More...
 
void reset (const label regioni, const labelUList &regions, const label patchID=-1, const bool syncCouples=true)
 Use the cells of cells corresponding to where region == regioni. More...
 
void setCellSubset (const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. Same as reset() More...
 
void setCellSubset (const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. Same as reset() More...
 
void setCellSubset (const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. Same as reset() More...
 
void setCellSubset (const label regioni, const labelUList &regions, const label patchID=-1, const bool syncPar=true)
 Use the specified subset of cells. Same as reset() More...
 
template<class Type >
tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &, const bool allowUnmapped=false) const
 Map volume internal (dimensioned) field Optional unmapped argument (currently unused) More...
 
template<class Type >
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const bool allowUnmapped=false) const
 Map volume field. More...
 
template<class Type >
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const bool allowUnmapped=false) const
 Map surface field. More...
 
template<class Type >
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate (const GeometricField< Type, pointPatchField, pointMesh > &, const bool allowUnmapped=false) const
 Map point field. More...
 

Static Public Member Functions

template<class Type >
static tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
 Map volume internal (dimensioned) field. More...
 
template<class Type >
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped=false)
 Map volume field. More...
 
template<class Type >
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap)
 Map surface field. More...
 
template<class Type >
static tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate (const GeometricField< Type, pointPatchField, pointMesh > &, const pointMesh &sMesh, const labelUList &patchMap, const labelUList &pointMap)
 Map point field. More...
 

Static Public Attributes

static word exposedPatchName
 Name for exposed internal faces (default: oldInternalFaces) More...
 

Protected Member Functions

bool checkHasSubMesh () const
 FatalError if subset has not been performed. More...
 
 fvMeshSubset (const fvMeshSubset &)=delete
 No copy construct. More...
 
void operator= (const fvMeshSubset &)=delete
 No copy assignment. More...
 

Detailed Description

Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the subMesh) with mapping lists for points, faces, and cells.

Can be constructed or reset to subset on the list of selected cells, which it creates as subMesh consisting only of the desired cells, with the mapping list for points, faces, and cells.

Places all exposed internal faces into either

  • a user supplied patch
  • a newly created patch "oldInternalFaces"
  • reset() does coupled patch subsetting as well. If it detects a face on a coupled patch 'losing' its neighbour it will move the face into the oldInternalFaces patch.
  • if a user supplied patch is used, it is up to the destination patchField to handle exposed internal faces (mapping from face -1). If not provided the default is to assign the internalField. All the basic patch field types (e.g. fixedValue) will give a warning and preferably derived patch field types should be used that know how to handle exposed faces (e.g. use uniformFixedValue instead of fixedValue)
Source files

Definition at line 75 of file fvMeshSubset.H.

Constructor & Destructor Documentation

◆ fvMeshSubset() [1/7]

fvMeshSubset ( const fvMeshSubset )
protecteddelete

No copy construct.

◆ fvMeshSubset() [2/7]

fvMeshSubset ( const fvMesh baseMesh)
explicit

Construct using the entire mesh (no subset)

Definition at line 394 of file fvMeshSubset.C.

◆ fvMeshSubset() [3/7]

fvMeshSubset ( const fvMesh baseMesh,
const Foam::zero   
)

Construct a zero-sized subset mesh, non-processor patches only.

Definition at line 406 of file fvMeshSubset.C.

References fvMeshSubset::reset().

Here is the call graph for this function:

◆ fvMeshSubset() [4/7]

fvMeshSubset ( const fvMesh baseMesh,
const bitSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

Construct for a cell-subset of the given mesh.

See reset() for more details.

Definition at line 415 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [5/7]

fvMeshSubset ( const fvMesh baseMesh,
const labelUList selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

Construct for a cell-subset of the given mesh.

See reset() for more details.

Definition at line 429 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [6/7]

fvMeshSubset ( const fvMesh baseMesh,
const labelHashSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

Construct for a cell-subset of the given mesh.

See reset() for more details.

Definition at line 443 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [7/7]

fvMeshSubset ( const fvMesh baseMesh,
const label  regioni,
const labelUList regions,
const label  patchID = -1,
const bool  syncPar = true 
)

Construct for a cell-subset of the given mesh.

See reset() for more details.

Definition at line 457 of file fvMeshSubset.C.

References patchID.

Member Function Documentation

◆ checkHasSubMesh()

bool checkHasSubMesh ( ) const
protected

FatalError if subset has not been performed.

Definition at line 106 of file fvMeshSubset.C.

References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and Foam::nl.

Here is the call graph for this function:

◆ operator=()

void operator= ( const fvMeshSubset )
protecteddelete

No copy assignment.

◆ baseMesh()

const Foam::fvMesh & baseMesh ( ) const
inlinenoexcept

Original mesh.

Definition at line 23 of file fvMeshSubsetI.H.

Referenced by zoneSubSet::baseMesh(), vtkWrite::writeVolFieldsImpl(), and ensightWrite::writeVolFieldsImpl().

Here is the caller graph for this function:

◆ mesh()

const Foam::fvMesh & mesh ( ) const
inlinenoexcept

Return baseMesh or subMesh, depending on the current state.

Definition at line 29 of file fvMeshSubsetI.H.

◆ hasSubMesh()

bool hasSubMesh ( ) const
inlinenoexcept

Have subMesh?

Definition at line 35 of file fvMeshSubsetI.H.

Referenced by fvMeshSubsetProxy::interpolate(), and fvMeshSubsetProxy::interpolateInternal().

Here is the caller graph for this function:

◆ subMesh() [1/2]

const Foam::fvMesh & subMesh ( ) const
inline

Return reference to subset mesh.

Definition at line 41 of file fvMeshSubsetI.H.

Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), marchingCells::initialise(), fvMeshSubsetProxy::mesh(), structuredRenumber::renumber(), and momentumError::write().

Here is the caller graph for this function:

◆ subMesh() [2/2]

Foam::fvMesh & subMesh ( )
inline

Return reference to subset mesh.

Definition at line 49 of file fvMeshSubsetI.H.

◆ pointMap()

const Foam::labelList & pointMap ( ) const
inline

Return point map.

Definition at line 57 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ faceMap()

const Foam::labelList & faceMap ( ) const
inline

Return face map.

Definition at line 65 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute(), and marchingCells::initialise().

Here is the caller graph for this function:

◆ faceFlipMap()

const Foam::labelList & faceFlipMap ( ) const
inline

Return face map with sign to encode flipped faces.

Definition at line 73 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ cellMap()

const Foam::labelList & cellMap ( ) const
inline

Return cell map.

Definition at line 84 of file fvMeshSubsetI.H.

Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), and structuredRenumber::renumber().

Here is the caller graph for this function:

◆ patchMap()

const Foam::labelList & patchMap ( ) const
inline

Return patch map.

Definition at line 92 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ clear()

void clear ( )

Reset subMesh and all maps.

Definition at line 473 of file fvMeshSubset.C.

◆ reset() [1/7]

void reset ( )

Reset subMesh and all maps. Same as clear()

Definition at line 485 of file fvMeshSubset.C.

References clear().

Referenced by fvMeshSubset::fvMeshSubset(), and fvMeshSubset::setCellSubset().

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

◆ reset() [2/7]

void reset ( const Foam::zero  )

Reset to a zero-sized subset mesh, non-processor patches only.

Definition at line 516 of file fvMeshSubset.C.

References clear(), forAll, Foam::identity(), polyBoundaryMesh::nNonProcessor(), IOobjectOption::NO_READ, and IOobjectOption::NO_WRITE.

Here is the call graph for this function:

◆ reset() [3/7]

void reset ( autoPtr< fvMesh > &&  subMeshPtr,
labelList &&  pointMap,
labelList &&  faceMap,
labelList &&  cellMap,
labelList &&  patchMap 
)

Reset from components.

Parameters
subMeshPtrMesh subset
pointMapPoint mapping
faceMapFace mapping
cellMapCell mapping
patchMapPatch mapping

Definition at line 492 of file fvMeshSubset.C.

References clear(), and Foam::faceMap().

Here is the call graph for this function:

◆ reset() [4/7]

void reset ( const bitSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

◆ reset() [5/7]

void reset ( const labelUList selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

Use the specified subset of cells.

Definition at line 1232 of file fvMeshSubset.C.

References Foam::BitSetOps::create(), and patchID.

Here is the call graph for this function:

◆ reset() [6/7]

void reset ( const labelHashSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)

Use the specified subset of cells.

Definition at line 1248 of file fvMeshSubset.C.

References Foam::BitSetOps::create(), and patchID.

Here is the call graph for this function:

◆ reset() [7/7]

void reset ( const label  regioni,
const labelUList regions,
const label  patchID = -1,
const bool  syncCouples = true 
)

Use the cells of cells corresponding to where region == regioni.

Definition at line 1264 of file fvMeshSubset.C.

References Foam::BitSetOps::create(), and patchID.

Here is the call graph for this function:

◆ setCellSubset() [1/4]

void setCellSubset ( const bitSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)
inline

Use the specified subset of cells. Same as reset()

Definition at line 374 of file fvMeshSubset.H.

References patchID, and fvMeshSubset::reset().

Here is the call graph for this function:

◆ setCellSubset() [2/4]

void setCellSubset ( const labelUList selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)
inline

Use the specified subset of cells. Same as reset()

Definition at line 387 of file fvMeshSubset.H.

References patchID, and fvMeshSubset::reset().

Here is the call graph for this function:

◆ setCellSubset() [3/4]

void setCellSubset ( const labelHashSet selectedCells,
const label  patchID = -1,
const bool  syncPar = true 
)
inline

Use the specified subset of cells. Same as reset()

Definition at line 400 of file fvMeshSubset.H.

References patchID, and fvMeshSubset::reset().

Here is the call graph for this function:

◆ setCellSubset() [4/4]

void setCellSubset ( const label  regioni,
const labelUList regions,
const label  patchID = -1,
const bool  syncPar = true 
)
inline

Use the specified subset of cells. Same as reset()

Definition at line 413 of file fvMeshSubset.H.

References patchID, and fvMeshSubset::reset().

Here is the call graph for this function:

◆ interpolate() [1/8]

static tmp<DimensionedField<Type, volMesh> > interpolate ( const DimensionedField< Type, volMesh > &  ,
const fvMesh sMesh,
const labelUList cellMap 
)
static

Map volume internal (dimensioned) field.

Referenced by momentumError::calcMomentError(), momentumError::divDevRhoReff(), and fvMeshSubsetProxy::interpolate().

Here is the caller graph for this function:

◆ interpolate() [2/8]

static tmp<GeometricField<Type, fvPatchField, volMesh> > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  ,
const fvMesh sMesh,
const labelUList patchMap,
const labelUList cellMap,
const labelUList faceMap,
const bool  allowUnmapped = false 
)
static

Map volume field.

Optionally allow unmapped faces not to produce a warning

◆ interpolate() [3/8]

static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate ( const GeometricField< Type, fvsPatchField, surfaceMesh > &  ,
const fvMesh sMesh,
const labelUList patchMap,
const labelUList cellMap,
const labelUList faceMap 
)
static

Map surface field.

Optionally negates value if flipping a face (from exposing an internal face)

◆ interpolate() [4/8]

static tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate ( const GeometricField< Type, pointPatchField, pointMesh > &  ,
const pointMesh sMesh,
const labelUList patchMap,
const labelUList pointMap 
)
static

Map point field.

◆ interpolate() [5/8]

tmp<DimensionedField<Type, volMesh> > interpolate ( const DimensionedField< Type, volMesh > &  ,
const bool  allowUnmapped = false 
) const

Map volume internal (dimensioned) field Optional unmapped argument (currently unused)

◆ interpolate() [6/8]

tmp<GeometricField<Type, fvPatchField, volMesh> > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  ,
const bool  allowUnmapped = false 
) const

Map volume field.

Optionally allow unmapped faces not to produce a warning

◆ interpolate() [7/8]

tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate ( const GeometricField< Type, fvsPatchField, surfaceMesh > &  ,
const bool  allowUnmapped = false 
) const

Map surface field.

Optionally allow unmapped faces not to produce a warning (currently unused)

◆ interpolate() [8/8]

tmp<GeometricField<Type, pointPatchField, pointMesh> > interpolate ( const GeometricField< Type, pointPatchField, pointMesh > &  ,
const bool  allowUnmapped = false 
) const

Map point field.

Optionally allow unmapped points not to produce a warning (currently unused)

Member Data Documentation

◆ exposedPatchName

Foam::word exposedPatchName
static

Name for exposed internal faces (default: oldInternalFaces)

Definition at line 164 of file fvMeshSubset.H.

Referenced by marchingCells::initialise().


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