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... | |
void | setCellSubset (const labelHashSet &globalCellMap, const label patchID=-1) |
Set the subset. Create "oldInternalFaces" patch for exposed. More... | |
void | setLargeCellSubset (const labelList ®ion, 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 ®ion, const label currentRegion, const bool syncCouples=true) const |
Two step subsetting. More... | |
void | setLargeCellSubset (const labelList ®ion, const label currentRegion, const labelList &exposedFaces, const labelList &patchIDs, const bool syncCouples=true) |
For every exposed face (from above getExposedFaces) More... | |
const fvMesh & | baseMesh () const |
Original mesh. More... | |
bool | hasSubMesh () const |
Have subMesh? More... | |
const fvMesh & | subMesh () const |
Return reference to subset mesh. More... | |
fvMesh & | subMesh () |
const labelList & | pointMap () const |
Return point map. More... | |
const labelList & | faceMap () const |
Return face map. More... | |
const labelList & | faceFlipMap () const |
Return face map with sign to encode flipped faces. More... | |
const labelList & | cellMap () const |
Return cell map. More... | |
const labelList & | patchMap () 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 field. More... | |
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
Definition at line 73 of file fvMeshSubset.H.
|
explicit |
Construct given a mesh to subset.
Definition at line 395 of file fvMeshSubset.C.
References fvMeshSubset::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(), fvMeshSubset::baseMesh(), polyMesh::boundaryMesh(), primitiveMesh::cells(), IOobject::clone(), Foam::endl(), polyMesh::faceNeighbour(), polyMesh::faceOwner(), polyMesh::faces(), Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, HashTable< T, Key, Hash >::found(), HashTable< T, Key, Hash >::insert(), HashTable< T, label, Hash< label > >::insert(), primitiveMesh::isInternalFace(), Foam::labelMax, Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, patchi, polyMesh::points(), Foam::Pout, fvMeshSubset::setLargeCellSubset(), List< T >::setSize(), List< T >::size(), UPtrList< T >::size(), HashTable< T, label, Hash< label > >::size(), HashTable< T, Key, Hash >::size(), Foam::sort(), timeName, HashTable< T, Key, Hash >::toc(), polyBoundaryMesh::whichPatch(), and Foam::xferMove().
Referenced by fvMeshSubset::fvMeshSubset().
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(), fvMeshSubset::baseMesh(), polyMesh::boundaryMesh(), primitiveMesh::cells(), IOobject::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(), primitiveMesh::nInternalFaces(), IOobject::NO_READ, IOobject::NO_WRITE, UPstream::nProcs(), UPstream::parRun(), patchi, patchNames(), polyMesh::points(), Foam::reduce(), Pstream::scatterList(), List< T >::setSize(), List< T >::size(), UPtrList< T >::size(), polyPatch::start(), timeName, and Foam::xferMove().
Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), fvMeshSubset::getExposedFaces(), structuredRenumber::renumber(), fvMeshSubset::setCellSubset(), and fvMeshSubset::setLargeCellSubset().
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 fvMeshSubset::baseMesh(), forAllConstIter(), fvMeshSubset::getExposedFaces(), and fvMeshSubset::setLargeCellSubset().
Foam::labelList getExposedFaces | ( | const labelList & | region, |
const label | currentRegion, | ||
const bool | syncCouples = true |
||
) | const |
Two step subsetting.
Get labels of exposed faces.
These are
Definition at line 1405 of file fvMeshSubset.C.
References fvMeshSubset::baseMesh(), removeCells::getExposedFaces(), and fvMeshSubset::setLargeCellSubset().
Referenced by fvMeshSubset::setLargeCellSubset().
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 fvMeshSubset::baseMesh(), Foam::identity(), polyTopoChange::makeMesh(), Foam::name(), IOobject::NO_READ, IOobject::NO_WRITE, removeCells::setRefinement(), and timeName.
|
inline |
Original mesh.
Definition at line 217 of file fvMeshSubset.H.
References fvMeshSubset::cellMap(), fvMeshSubset::faceFlipMap(), fvMeshSubset::faceMap(), fvMeshSubset::hasSubMesh(), fvMeshSubset::interpolate(), fvMeshSubset::patchMap(), fvMeshSubset::pointMap(), and fvMeshSubset::subMesh().
Referenced by fvMeshSubset::getExposedFaces(), fvMeshDistribute::printFieldInfo(), fvMeshSubset::setCellSubset(), and fvMeshSubset::setLargeCellSubset().
bool hasSubMesh | ( | ) | const |
Have subMesh?
Definition at line 1466 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh().
const fvMesh & subMesh | ( | ) | const |
Return reference to subset mesh.
Definition at line 1472 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh(), structuredDecomp::decompose(), fvMeshDistribute::distribute(), fvMeshSubset::faceFlipMap(), vtkMesh::mesh(), and structuredRenumber::renumber().
fvMesh & subMesh | ( | ) |
Definition at line 1480 of file fvMeshSubset.C.
const labelList & pointMap | ( | ) | const |
Return point map.
Definition at line 1488 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh(), and fvMeshDistribute::distribute().
const labelList & faceMap | ( | ) | const |
Return face map.
Definition at line 1496 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh(), fvMeshDistribute::distribute(), and fvMeshSubset::faceFlipMap().
const labelList & faceFlipMap | ( | ) | const |
Return face map with sign to encode flipped faces.
Definition at line 1504 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().
const labelList & cellMap | ( | ) | const |
Return cell map.
Definition at line 1542 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh(), structuredDecomp::decompose(), fvMeshDistribute::distribute(), fvMeshSubset::faceFlipMap(), and structuredRenumber::renumber().
const labelList & patchMap | ( | ) | const |
Return patch map.
Definition at line 1550 of file fvMeshSubset.C.
Referenced by fvMeshSubset::baseMesh(), and fvMeshDistribute::distribute().
|
static |
Map volume field.
Definition at line 43 of file fvMeshSubsetInterpolate.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), DimensionedField< Type, GeoMesh >::dimensions(), forAll, DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), IOobject::NO_READ, IOobject::NO_WRITE, DimensionedField< Type, GeoMesh >::null(), patchi, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), GeometricField< Type, PatchField, GeoMesh >::ref(), List< T >::size(), fvPatch::size(), fvPatch::start(), fvMesh::time(), and Time::timeName().
Referenced by fvMeshSubset::baseMesh(), vtkMesh::interpolate(), fvMeshSubset::interpolate(), and fvMeshDistribute::printFieldInfo().
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) | const |
Definition at line 161 of file fvMeshSubsetInterpolate.C.
References Foam::faceMap(), Foam::interpolate(), and fvMeshSubset::interpolate().
|
static |
Map surface field. Optionally negates value if flipping.
a face (from exposing an internal face)
Definition at line 178 of file fvMeshSubsetInterpolate.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::faceMap(), forAll, GeometricField< Type, PatchField, GeoMesh >::internalField(), fvMeshSubset::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), primitiveMesh::nInternalFaces(), IOobject::NO_READ, IOobject::NO_WRITE, DimensionedField< Type, GeoMesh >::null(), fvPatch::patch(), patchi, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), List< T >::size(), fvPatch::size(), fvPatch::start(), fvMesh::time(), Time::timeName(), and polyPatch::whichFace().
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 fvMeshSubset::interpolate().
|
static |
Map point field.
Definition at line 364 of file fvMeshSubsetInterpolate.C.
References pointMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), DimensionedField< Type, GeoMesh >::dimensions(), HashTableCore::end(), HashTable< T, Key, Hash >::find(), forAll, fvMeshSubset::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), pointPatch::meshPoints(), IOobject::name(), IOobject::NO_READ, IOobject::NO_WRITE, DimensionedField< Type, GeoMesh >::null(), patchi, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), GeometricField< Type, PatchField, GeoMesh >::ref(), List< T >::size(), pointPatch::size(), pointMesh::thisDb(), IOobject::time(), and Time::timeName().
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | sf | ) | const |
Definition at line 494 of file fvMeshSubsetInterpolate.C.
References Foam::interpolate(), fvMeshSubset::interpolate(), and MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New().
|
static |
Map dimensioned field.
Definition at line 510 of file fvMeshSubsetInterpolate.C.
References DimensionedField< Type, GeoMesh >::dimensions(), fvMeshSubset::interpolate(), IOobject::name(), IOobject::NO_READ, IOobject::NO_WRITE, fvMesh::time(), and Time::timeName().
tmp< DimensionedField< Type, volMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | df | ) | const |
Definition at line 541 of file fvMeshSubsetInterpolate.C.
References Foam::interpolate().