fvMeshSubset Class Reference

Post-processing mesh subset tool. 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...

Public Member Functions

 fvMeshSubset (const fvMesh &)
 Construct given a mesh to subset. More...
 
 fvMeshSubset (const fvMeshSubset &)=delete
 Disallow default bitwise copy construction. More...
 
void setCellSubset (const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
 Set the subset. Create "oldInternalFaces" patch for exposed. More...
 
void setLargeCellSubset (const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
 Set the subset from all cells with region == currentRegion. More...
 
void setLargeCellSubset (const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
 setLargeCellSubset but with labelHashSet. More...
 
labelList getExposedFaces (const labelList &region, const label currentRegion, const bool syncCouples=true) const
 Two step subsetting. More...
 
void setLargeCellSubset (const labelList &region, const label currentRegion, const labelList &exposedFaces, const labelList &patchIDs, const bool syncCouples=true)
 For every exposed face (from above getExposedFaces) More...
 
const fvMeshbaseMesh () const
 Original mesh. More...
 
bool hasSubMesh () const
 Have subMesh? More...
 
const fvMeshsubMesh () const
 Return reference to subset mesh. More...
 
fvMeshsubMesh ()
 
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...
 
template<class Type >
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &) const
 
template<class Type >
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const bool negateIfFlipped=true) const
 
template<class Type >
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate (const GeometricField< Type, pointPatchField, pointMesh > &) const
 
template<class Type >
tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &) const
 
void operator= (const fvMeshSubset &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Member Functions

template<class Type >
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
 Map volume field. More...
 
template<class Type >
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap, const bool negateIfFlipped=true)
 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 labelList &patchMap, const labelList &pointMap)
 Map point field. More...
 
template<class Type >
static tmp< DimensionedField< Type, volMesh > > interpolate (const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelList &cellMap)
 Map dimensioned field. More...
 

Detailed Description

Post-processing mesh subset tool. 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 is for small subsets. Uses Maps to minimize memory.
  • setLargeCellSubset is for largish subsets (>10% of mesh). Uses labelLists instead.
  • setLargeCellSubset 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/2]

fvMeshSubset ( const fvMesh baseMesh)
explicit

Construct given a mesh to subset.

Definition at line 468 of file fvMeshSubset.C.

References fvMeshSubset::setCellSubset().

Here is the call graph for this function:

◆ fvMeshSubset() [2/2]

fvMeshSubset ( const fvMeshSubset )
delete

Disallow default bitwise copy construction.

Member Function Documentation

◆ setCellSubset()

◆ setLargeCellSubset() [1/3]

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

◆ setLargeCellSubset() [2/3]

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

setLargeCellSubset but with labelHashSet.

Definition at line 1479 of file fvMeshSubset.C.

References fvMeshSubset::baseMesh(), forAllConstIter(), fvMeshSubset::getExposedFaces(), and fvMeshSubset::setLargeCellSubset().

Here is the call graph for this function:

◆ getExposedFaces()

Foam::labelList getExposedFaces ( const labelList region,
const label  currentRegion,
const bool  syncCouples = true 
) const

Two step subsetting.

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 1496 of file fvMeshSubset.C.

References fvMeshSubset::baseMesh(), removeCells::getExposedFaces(), and fvMeshSubset::setLargeCellSubset().

Referenced by fvMeshSubset::setLargeCellSubset().

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

◆ setLargeCellSubset() [3/3]

void setLargeCellSubset ( const labelList region,
const label  currentRegion,
const labelList exposedFaces,
const labelList patchIDs,
const bool  syncCouples = true 
)

For every exposed face (from above getExposedFaces)

used supplied (existing!) patch

Definition at line 1510 of file fvMeshSubset.C.

References fvMeshSubset::baseMesh(), Foam::identity(), polyTopoChange::makeMesh(), Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, removeCells::setRefinement(), and timeName.

Here is the call graph for this function:

◆ baseMesh()

◆ hasSubMesh()

bool hasSubMesh ( ) const

Have subMesh?

Definition at line 1557 of file fvMeshSubset.C.

Referenced by fvMeshSubset::baseMesh().

Here is the caller graph for this function:

◆ subMesh() [1/2]

◆ subMesh() [2/2]

fvMesh & subMesh ( )

Definition at line 1571 of file fvMeshSubset.C.

◆ pointMap()

const labelList & pointMap ( ) const

Return point map.

Definition at line 1579 of file fvMeshSubset.C.

Referenced by fvMeshSubset::baseMesh(), and fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ faceMap()

const labelList & faceMap ( ) const

◆ faceFlipMap()

const labelList & faceFlipMap ( ) const

Return face map with sign to encode flipped faces.

Definition at line 1595 of file fvMeshSubset.C.

References fvMeshSubset::cellMap(), fvMeshSubset::faceMap(), polyMesh::faceOwner(), primitiveMesh::nInternalFaces(), List< T >::size(), and fvMeshSubset::subMesh().

Referenced by fvMeshSubset::baseMesh(), and fvMeshDistribute::distribute().

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

◆ cellMap()

const labelList & cellMap ( ) const

Return cell map.

Definition at line 1633 of file fvMeshSubset.C.

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

Here is the caller graph for this function:

◆ patchMap()

const labelList & patchMap ( ) const

Return patch map.

Definition at line 1641 of file fvMeshSubset.C.

Referenced by waveAlphaFvPatchScalarField::alphan(), fvMeshSubset::baseMesh(), fvMeshDistribute::distribute(), wavePressureFvPatchScalarField::pn(), and waveVelocityFvPatchVectorField::Un().

Here is the caller graph for this function:

◆ interpolate() [1/8]

◆ interpolate() [2/8]

tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const

Definition at line 162 of file fvMeshSubsetInterpolate.C.

References Foam::faceMap(), Foam::interpolate(), and fvMeshSubset::interpolate().

Here is the call graph for this function:

◆ interpolate() [3/8]

◆ interpolate() [4/8]

tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate ( const GeometricField< Type, fvsPatchField, surfaceMesh > &  sf,
const bool  negateIfFlipped = true 
) const

Definition at line 346 of file fvMeshSubsetInterpolate.C.

References Foam::faceMap(), Foam::interpolate(), and fvMeshSubset::interpolate().

Here is the call graph for this function:

◆ interpolate() [5/8]

◆ interpolate() [6/8]

tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate ( const GeometricField< Type, pointPatchField, pointMesh > &  sf) const

Definition at line 497 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate(), fvMeshSubset::interpolate(), and MeshObject< polyMesh, PatchMeshObject, pointMesh >::New().

Here is the call graph for this function:

◆ interpolate() [7/8]

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

◆ interpolate() [8/8]

tmp< DimensionedField< Type, volMesh > > interpolate ( const DimensionedField< Type, volMesh > &  df) const

Definition at line 544 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate().

Here is the call graph for this function:

◆ operator=()

void operator= ( const fvMeshSubset )
delete

Disallow default bitwise assignment.

Referenced by fvMeshSubset::baseMesh().

Here is the caller graph for this function:

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