fvMeshSubset Class Reference

Given the original mesh and the list of selected cells, it creates the mesh consisting only of the desired cells, with the mapping list for points, faces, and cells. More...

Collaboration diagram for fvMeshSubset:
[legend]

Public Member Functions

 fvMeshSubset (const fvMesh &baseMesh)
 Construct given a mesh to subset. 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 maps and subsetting. More...
 
void setCellSubset (const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 Define cell subset based on the selectedCells. More...
 
void setCellSubset (const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
 
void setCellSubset (const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
 
void setCellSubset (const label regioni, const labelUList &regions, const label patchID=-1, const bool syncCouples=true)
 Define cell subset, using the cells for which region == regioni. More...
 
labelList getExposedFaces (const bitSet &selectedCells, const bool syncCouples=true) const
 Get labels of exposed faces. More...
 
labelList getExposedFaces (const label regioni, const labelUList &regions, const bool syncCouples=true) const
 Get labels of exposed faces. More...
 
void setCellSubset (const bitSet &selectedCells, const labelList &exposedFaces, const labelList &patchIDs, const bool syncCouples=true)
 For every exposed face (from above getExposedFaces) More...
 
void setCellSubset (const label regioni, const labelList &regions, const labelList &exposedFaces, const labelList &patchIDs, const bool syncCouples=true)
 For every exposed face (from above getExposedFaces) More...
 
template<class Type >
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const bool allowUnmapped=false) const
 
template<class Type >
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. More...
 
template<class Type >
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. More...
 
template<class Type >
tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &, const bool allowUnmapped=false) const
 Map Dimensioned. Optional unmapped argument not used. More...
 
void setLargeCellSubset (const labelUList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
 Deprecated(2018-07) old method name and old parameter ordering. More...
 
void setLargeCellSubset (const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
 Deprecated(2018-07) old method name. More...
 
void setLargeCellSubset (const labelList &regions, const label regioni, const labelList &exposedFaces, const labelList &patchIDs, const bool syncCouples=true)
 Deprecated(2018-07) method. More...
 

Static Public Member Functions

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. Optionally allow unmapped faces not to produce. 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. Optionally negates value if flipping. 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...
 
template<class Type >
static tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
 Map dimensioned field. More...
 

Static Public Attributes

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

Detailed Description

Given the original mesh and the list of selected cells, it creates the mesh consisting only of the desired cells, with the mapping list for points, faces, and cells.

Puts all exposed internal faces into either

  • a user supplied patch
  • a newly created patch "oldInternalFaces"
  • setCellSubset 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 73 of file fvMeshSubset.H.

Constructor & Destructor Documentation

◆ fvMeshSubset() [1/5]

fvMeshSubset ( const fvMesh baseMesh)
explicit

Construct given a mesh to subset.

Definition at line 474 of file fvMeshSubset.C.

◆ fvMeshSubset() [2/5]

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 setCellSubset() for more details.

Definition at line 487 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [3/5]

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 setCellSubset() for more details.

Definition at line 515 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [4/5]

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 setCellSubset() for more details.

Definition at line 501 of file fvMeshSubset.C.

References patchID.

◆ fvMeshSubset() [5/5]

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 setCellSubset() for more details.

Definition at line 529 of file fvMeshSubset.C.

References patchID.

Member Function Documentation

◆ baseMesh()

const Foam::fvMesh & baseMesh ( ) const
inlinenoexcept

Original mesh.

Definition at line 30 of file fvMeshSubsetI.H.

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

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 36 of file fvMeshSubsetI.H.

◆ hasSubMesh()

bool hasSubMesh ( ) const
inlinenoexcept

Have subMesh?

Definition at line 42 of file fvMeshSubsetI.H.

References bool.

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 48 of file fvMeshSubsetI.H.

Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), 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 56 of file fvMeshSubsetI.H.

◆ pointMap()

const Foam::labelList & pointMap ( ) const
inline

Return point map.

Definition at line 64 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 72 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute().

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 80 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 91 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 99 of file fvMeshSubsetI.H.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ clear()

void clear ( )

Reset maps and subsetting.

Definition at line 545 of file fvMeshSubset.C.

◆ setCellSubset() [1/6]

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

Define cell subset based on the selectedCells.

Create "oldInternalFaces" patch for exposed internal faces (patchID==-1) or use supplied patch. Handles coupled patches if necessary by making coupled patch face part of patchID (so uncoupled)

Definition at line 558 of file fvMeshSubset.C.

References Foam::abort(), polyMesh::boundaryMesh(), primitiveMesh::cells(), clear(), PtrList< T >::clone(), Foam::endl(), polyMesh::faceNeighbour(), polyMesh::faceOwner(), polyMesh::faces(), Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, Pstream::gatherList(), Foam::identity(), primitiveMesh::isInternalFace(), Foam::labelMax, Pstream::listCombineGather(), Pstream::listCombineScatter(), UPstream::myProcNo(), Foam::name(), polyBoundaryMesh::names(), autoPtr< T >::New(), primitiveMesh::nInternalFaces(), IOobject::NO_READ, IOobject::NO_WRITE, UPstream::nProcs(), UPstream::parRun(), patchID, patchNames(), polyMesh::points(), Foam::reduce(), List< T >::resize(), Pstream::scatterList(), List< T >::setSize(), bitSet::sortedToc(), polyPatch::start(), bitSet::test(), timeName, and Foam::Zero.

Referenced by fvMeshSubset::setLargeCellSubset().

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

◆ setCellSubset() [2/6]

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

Define cell subset, using the specified cells to define the selection

Definition at line 1196 of file fvMeshSubset.C.

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

Here is the call graph for this function:

◆ setCellSubset() [3/6]

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

Define cell subset, using the specified cells labelHashSet to define the selection

Definition at line 1212 of file fvMeshSubset.C.

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

Here is the call graph for this function:

◆ setCellSubset() [4/6]

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

Define cell subset, using the cells for which region == regioni.

Definition at line 1228 of file fvMeshSubset.C.

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

Here is the call graph for this function:

◆ getExposedFaces() [1/2]

Foam::labelList getExposedFaces ( const bitSet selectedCells,
const bool  syncCouples = true 
) const

Get labels of exposed faces.

These are

  • internal faces that become boundary faces
  • coupled faces that become uncoupled (since one of the sides gets deleted)

Definition at line 1245 of file fvMeshSubset.C.

References removeCells::getExposedFaces().

Here is the call graph for this function:

◆ getExposedFaces() [2/2]

Foam::labelList getExposedFaces ( const label  regioni,
const labelUList regions,
const bool  syncCouples = true 
) const

Get labels of exposed faces.

These are

  • internal faces that become boundary faces
  • coupled faces that become uncoupled (since one of the sides gets deleted)

Definition at line 1257 of file fvMeshSubset.C.

References removeCells::getExposedFaces().

Here is the call graph for this function:

◆ setCellSubset() [5/6]

void setCellSubset ( const bitSet selectedCells,
const labelList exposedFaces,
const labelList patchIDs,
const bool  syncCouples = true 
)

For every exposed face (from above getExposedFaces)

used supplied (existing!) patch

Definition at line 1270 of file fvMeshSubset.C.

◆ setCellSubset() [6/6]

void setCellSubset ( const label  regioni,
const labelList regions,
const labelList exposedFaces,
const labelList patchIDs,
const bool  syncCouples = true 
)

For every exposed face (from above getExposedFaces)

used supplied (existing!) patch

Definition at line 1288 of file fvMeshSubset.C.

◆ interpolate() [1/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

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

Here is the caller graph for this function:

◆ interpolate() [2/8]

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

◆ 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]

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 (not currently used)

◆ interpolate() [5/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() [6/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 (not currently used)

◆ interpolate() [7/8]

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

Map dimensioned field.

◆ interpolate() [8/8]

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

Map Dimensioned. Optional unmapped argument not used.

◆ setLargeCellSubset() [1/3]

void setLargeCellSubset ( const labelUList region,
const label  currentRegion,
const label  patchID = -1,
const bool  syncCouples = true 
)
inline

Deprecated(2018-07) old method name and old parameter ordering.

Deprecated:
(2018-07) - use setCellSubset() method

Definition at line 432 of file fvMeshSubset.H.

References Foam::Info, Foam::nl, patchID, and fvMeshSubset::setCellSubset().

Here is the call graph for this function:

◆ setLargeCellSubset() [2/3]

void setLargeCellSubset ( const labelHashSet globalCellMap,
const label  patchID = -1,
const bool  syncPar = true 
)
inline

Deprecated(2018-07) old method name.

Deprecated:
(2018-07) - use setCellSubset() method

Definition at line 455 of file fvMeshSubset.H.

References Foam::Info, Foam::nl, patchID, and fvMeshSubset::setCellSubset().

Here is the call graph for this function:

◆ setLargeCellSubset() [3/3]

void setLargeCellSubset ( const labelList regions,
const label  regioni,
const labelList exposedFaces,
const labelList patchIDs,
const bool  syncCouples = true 
)
inline

Deprecated(2018-07) method.

For every exposed face (from getExposedFaces) use supplied (existing!) patch ids

Deprecated:
(2018-07) - use setCellSubset() method

Definition at line 473 of file fvMeshSubset.H.

References Foam::Info, Foam::nl, and fvMeshSubset::setCellSubset().

Here is the call graph for this function:

Member Data Documentation

◆ exposedPatchName

Foam::word exposedPatchName
static

Name for exposed internal faces (default: oldInternalFaces)

Definition at line 161 of file fvMeshSubset.H.


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