Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
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...

Collaboration diagram for fvMeshSubset:
Collaboration graph
[legend]

Public Member Functions

 fvMeshSubset (const fvMesh &)
 Construct given a mesh to subset. More...
 
void setCellSubset (const labelHashSet &globalCellMap, const label patchID=-1)
 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
 

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 fields. More...
 

Private Member Functions

bool checkCellSubset () const
 Check if subset has been performed. More...
 
void doCoupledPatches (const bool syncPar, labelList &nCellsUsingFace) const
 Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'. More...
 
void subsetZones ()
 Create zones for submesh. More...
 
labelList getCellsToRemove (const labelList &region, const label currentRegion) const
 Helper: extract cells-to-remove from cells-to-keep. More...
 
 fvMeshSubset (const fvMeshSubset &)
 Disallow default bitwise copy construct. More...
 
void operator= (const fvMeshSubset &)
 Disallow default bitwise assignment. More...
 

Static Private Member Functions

static void markPoints (const labelList &, Map< label > &)
 Mark points in Map. More...
 
static void markPoints (const labelList &, labelList &)
 Mark points (with 0) in labelList. More...
 
static labelList subset (const label nElems, const labelList &selectedElements, const labelList &subsetMap)
 Subset of subset. More...
 

Private Attributes

const fvMeshbaseMesh_
 Mesh to subset from. More...
 
autoPtr< fvMeshfvMeshSubsetPtr_
 Subset mesh pointer. More...
 
labelList pointMap_
 Point mapping array. More...
 
labelList faceMap_
 Face mapping array. More...
 
labelList cellMap_
 Cell mapping array. More...
 
labelList patchMap_
 Patch mapping array. More...
 
autoPtr< labelListfaceFlipMapPtr_
 

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

Source files

Definition at line 73 of file fvMeshSubset.H.

Constructor & Destructor Documentation

◆ fvMeshSubset() [1/2]

fvMeshSubset ( const fvMeshSubset )
private

Disallow default bitwise copy construct.

◆ fvMeshSubset() [2/2]

fvMeshSubset ( const fvMesh baseMesh)
explicit

Construct given a mesh to subset.

Definition at line 395 of file fvMeshSubset.C.

Member Function Documentation

◆ checkCellSubset()

bool checkCellSubset ( ) const
private

Check if subset has been performed.

Definition at line 48 of file fvMeshSubset.C.

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

Here is the call graph for this function:

◆ markPoints() [1/2]

void markPoints ( const labelList curPoints,
Map< label > &  pointMap 
)
staticprivate

Mark points in Map.

Definition at line 67 of file fvMeshSubset.C.

References forAll.

◆ markPoints() [2/2]

void markPoints ( const labelList curPoints,
labelList pointMap 
)
staticprivate

Mark points (with 0) in labelList.

Definition at line 81 of file fvMeshSubset.C.

References forAll.

◆ doCoupledPatches()

void doCoupledPatches ( const bool  syncPar,
labelList nCellsUsingFace 
) const
private

Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.

Definition at line 96 of file fvMeshSubset.C.

References Foam::endl(), PstreamBuffers::finishedSends(), forAll, Foam::Info, processorPolyPatch::neighbProcNo(), UPstream::nonBlocking, UPstream::parRun(), Foam::reduce(), polyPatch::start(), and cyclicPolyPatch::transformGlobalFace().

Here is the call graph for this function:

◆ subset()

labelList subset ( const label  nElems,
const labelList selectedElements,
const labelList subsetMap 
)
staticprivate

Subset of subset.

Definition at line 204 of file fvMeshSubset.C.

References forAll, and n.

◆ subsetZones()

void subsetZones ( )
private

Create zones for submesh.

Definition at line 243 of file fvMeshSubset.C.

References Foam::faceMap(), faceZone::flipMap(), forAll, zone::name(), nPoints, and Foam::subset().

Here is the call graph for this function:

◆ getCellsToRemove()

Foam::labelList getCellsToRemove ( const labelList region,
const label  currentRegion 
) const
private

Helper: extract cells-to-remove from cells-to-keep.

Definition at line 360 of file fvMeshSubset.C.

References forAll.

◆ operator=()

void operator= ( const fvMeshSubset )
private

Disallow default bitwise assignment.

◆ setCellSubset()

void setCellSubset ( const labelHashSet globalCellMap,
const label  patchID = -1 
)

Set the subset. Create "oldInternalFaces" patch for exposed.

internal faces (patchID==-1) or use supplied patch. Does not handle coupled patches correctly if only one side gets deleted.

Definition at line 410 of file fvMeshSubset.C.

References Foam::abort(), IOobject::clone(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, HashTable::found(), Foam::labelMax, Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, Foam::Pout, List::setSize(), List::size(), PtrList::size(), Foam::sort(), timeName, HashTable::toc(), polyBoundaryMesh::whichPatch(), and Foam::xferMove().

Here is the call graph for this function:

◆ setLargeCellSubset() [1/3]

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.

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

Definition at line 795 of file fvMeshSubset.C.

References Foam::abort(), IOobject::clone(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, Pstream::gatherList(), Foam::identity(), Foam::labelMax, Pstream::listCombineGather(), Pstream::listCombineScatter(), UPstream::myProcNo(), Foam::name(), polyBoundaryMesh::names(), IOobject::NO_READ, IOobject::NO_WRITE, UPstream::nProcs(), UPstream::parRun(), patchNames(), Foam::reduce(), Pstream::scatterList(), List::setSize(), List::size(), PtrList::size(), polyPatch::start(), timeName, and Foam::xferMove().

Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), main(), vtkMesh::readUpdate(), regionRenumber(), and structuredRenumber::renumber().

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

◆ setLargeCellSubset() [2/3]

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

setLargeCellSubset but with labelHashSet.

Definition at line 1388 of file fvMeshSubset.C.

References forAllConstIter().

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

References removeCells::getExposedFaces().

Referenced by main().

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

References 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()

const fvMesh& baseMesh ( ) const
inline

Original mesh.

Definition at line 216 of file fvMeshSubset.H.

References fvMeshSubset::baseMesh_.

Referenced by fvMeshDistribute::sendFields(), subsetDimensionedFields(), subsetPointFields(), subsetSurfaceFields(), and subsetVolFields().

Here is the caller graph for this function:

◆ hasSubMesh()

bool hasSubMesh ( ) const

Have subMesh?

Definition at line 1466 of file fvMeshSubset.C.

◆ subMesh() [1/2]

fvMesh & subMesh ( ) const

Return reference to subset mesh.

Definition at line 1472 of file fvMeshSubset.C.

Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), main(), vtkMesh::mesh(), regionRenumber(), and structuredRenumber::renumber().

Here is the caller graph for this function:

◆ subMesh() [2/2]

fvMesh& subMesh ( )

◆ pointMap()

const labelList & pointMap ( ) const

Return point map.

Definition at line 1488 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ faceMap()

const labelList & faceMap ( ) const

Return face map.

Definition at line 1496 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ faceFlipMap()

const labelList & faceFlipMap ( ) const

Return face map with sign to encode flipped faces.

Definition at line 1504 of file fvMeshSubset.C.

References Foam::faceMap(), and List::size().

Referenced by 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 1542 of file fvMeshSubset.C.

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

Here is the caller graph for this function:

◆ patchMap()

const labelList & patchMap ( ) const

Return patch map.

Definition at line 1550 of file fvMeshSubset.C.

Referenced by fvMeshDistribute::distribute().

Here is the caller graph for this function:

◆ interpolate() [1/8]

tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  vf,
const fvMesh sMesh,
const labelList patchMap,
const labelList cellMap,
const labelList faceMap 
)
static

◆ interpolate() [2/8]

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

Definition at line 161 of file fvMeshSubsetInterpolate.C.

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

Here is the call graph for this function:

◆ interpolate() [3/8]

tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate ( const GeometricField< Type, fvsPatchField, surfaceMesh > &  vf,
const fvMesh sMesh,
const labelList patchMap,
const labelList cellMap,
const labelList faceMap,
const bool  negateIfFlipped = true 
)
static

◆ interpolate() [4/8]

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

Definition at line 344 of file fvMeshSubsetInterpolate.C.

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

Here is the call graph for this function:

◆ interpolate() [5/8]

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

◆ interpolate() [6/8]

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

Definition at line 494 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate(), MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), and sf().

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

Map dimensioned fields.

Definition at line 510 of file fvMeshSubsetInterpolate.C.

References DimensionedField::dimensions(), IOobject::NO_READ, IOobject::NO_WRITE, fvMesh::time(), and Time::timeName().

Here is the call graph for this function:

◆ interpolate() [8/8]

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

Definition at line 541 of file fvMeshSubsetInterpolate.C.

References Foam::interpolate().

Here is the call graph for this function:

Field Documentation

◆ baseMesh_

const fvMesh& baseMesh_
private

Mesh to subset from.

Definition at line 81 of file fvMeshSubset.H.

Referenced by fvMeshSubset::baseMesh().

◆ fvMeshSubsetPtr_

autoPtr<fvMesh> fvMeshSubsetPtr_
private

Subset mesh pointer.

Definition at line 84 of file fvMeshSubset.H.

◆ pointMap_

labelList pointMap_
private

Point mapping array.

Definition at line 87 of file fvMeshSubset.H.

◆ faceMap_

labelList faceMap_
private

Face mapping array.

Definition at line 90 of file fvMeshSubset.H.

◆ cellMap_

labelList cellMap_
private

Cell mapping array.

Definition at line 93 of file fvMeshSubset.H.

◆ patchMap_

labelList patchMap_
private

Patch mapping array.

Definition at line 96 of file fvMeshSubset.H.

◆ faceFlipMapPtr_

autoPtr<labelList> faceFlipMapPtr_
mutableprivate

Definition at line 98 of file fvMeshSubset.H.


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