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...
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 ®ions, const label patchID=-1, const bool syncPar=true) | |
Construct for a cell-subset of the given mesh. More... | |
const fvMesh & | baseMesh () const noexcept |
Original mesh. More... | |
const fvMesh & | mesh () const noexcept |
Return baseMesh or subMesh, depending on the current state. More... | |
bool | hasSubMesh () const noexcept |
Have subMesh? More... | |
const fvMesh & | subMesh () const |
Return reference to subset mesh. More... | |
fvMesh & | subMesh () |
Return reference to subset mesh. More... | |
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... | |
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 ®ions, 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 ®ions, 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... | |
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
Definition at line 75 of file fvMeshSubset.H.
|
protecteddelete |
No copy construct.
|
explicit |
Construct using the entire mesh (no subset)
Definition at line 394 of file fvMeshSubset.C.
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().
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 | ( | 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 | ( | 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 | ( | 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.
|
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.
|
protecteddelete |
No copy assignment.
|
inlinenoexcept |
Original mesh.
Definition at line 23 of file fvMeshSubsetI.H.
Referenced by zoneSubSet::baseMesh(), vtkWrite::writeVolFieldsImpl(), and ensightWrite::writeVolFieldsImpl().
|
inlinenoexcept |
Return baseMesh or subMesh, depending on the current state.
Definition at line 29 of file fvMeshSubsetI.H.
|
inlinenoexcept |
Have subMesh?
Definition at line 35 of file fvMeshSubsetI.H.
Referenced by fvMeshSubsetProxy::interpolate(), and fvMeshSubsetProxy::interpolateInternal().
|
inline |
Return reference to subset mesh.
Definition at line 41 of file fvMeshSubsetI.H.
Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), fvMeshSubsetProxy::mesh(), structuredRenumber::renumber(), and momentumError::write().
|
inline |
Return reference to subset mesh.
Definition at line 49 of file fvMeshSubsetI.H.
|
inline |
Return point map.
Definition at line 57 of file fvMeshSubsetI.H.
Referenced by fvMeshDistribute::distribute().
|
inline |
Return face map.
Definition at line 65 of file fvMeshSubsetI.H.
Referenced by fvMeshDistribute::distribute().
|
inline |
Return face map with sign to encode flipped faces.
Definition at line 73 of file fvMeshSubsetI.H.
Referenced by fvMeshDistribute::distribute().
|
inline |
Return cell map.
Definition at line 84 of file fvMeshSubsetI.H.
Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), and structuredRenumber::renumber().
|
inline |
Return patch map.
Definition at line 92 of file fvMeshSubsetI.H.
Referenced by fvMeshDistribute::distribute().
void clear | ( | ) |
Reset subMesh and all maps.
Definition at line 473 of file fvMeshSubset.C.
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().
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.
void reset | ( | autoPtr< fvMesh > && | subMeshPtr, |
labelList && | pointMap, | ||
labelList && | faceMap, | ||
labelList && | cellMap, | ||
labelList && | patchMap | ||
) |
Reset from components.
subMeshPtr | Mesh subset |
pointMap | Point mapping |
faceMap | Face mapping |
cellMap | Cell mapping |
patchMap | Patch mapping |
Definition at line 492 of file fvMeshSubset.C.
References clear(), and Foam::faceMap().
void reset | ( | const bitSet & | selectedCells, |
const label | patchID = -1 , |
||
const bool | syncPar = true |
||
) |
Use the specified subset of cells.
Definition at line 575 of file fvMeshSubset.C.
References Foam::abort(), Pstream::allGatherList(), UList< T >::begin(), polyMesh::boundaryMesh(), primitiveMesh::cells(), clear(), Foam::endl(), polyMesh::faceNeighbour(), polyMesh::faceOwner(), polyMesh::faces(), Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, Foam::identity(), primitiveMesh::isInternalFace(), Foam::labelMax, Pstream::listCombineReduce(), UPstream::myProcNo(), polyBoundaryMesh::names(), autoPtr< T >::New(), primitiveMesh::nInternalFaces(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, UPstream::nProcs(), UPstream::parRun(), patchID, patchNames(), polyMesh::points(), Foam::reduce(), List< T >::resize(), autoPtr< T >::set(), List< T >::setSize(), UList< T >::size(), UPtrList< T >::size(), bitSet::sortedToc(), polyPatch::start(), bitSet::test(), and Foam::Zero.
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.
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.
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.
|
inline |
Use the specified subset of cells. Same as reset()
Definition at line 374 of file fvMeshSubset.H.
References patchID, and fvMeshSubset::reset().
|
inline |
Use the specified subset of cells. Same as reset()
Definition at line 387 of file fvMeshSubset.H.
References patchID, and fvMeshSubset::reset().
|
inline |
Use the specified subset of cells. Same as reset()
Definition at line 400 of file fvMeshSubset.H.
References patchID, and fvMeshSubset::reset().
|
inline |
Use the specified subset of cells. Same as reset()
Definition at line 413 of file fvMeshSubset.H.
References patchID, and fvMeshSubset::reset().
|
static |
Map volume internal (dimensioned) field.
Referenced by momentumError::calcMomentError(), momentumError::divDevRhoReff(), and fvMeshSubsetProxy::interpolate().
|
static |
Map volume field.
Optionally allow unmapped faces not to produce a warning
|
static |
Map surface field.
Optionally negates value if flipping a face (from exposing an internal face)
|
static |
Map point field.
tmp<DimensionedField<Type, volMesh> > interpolate | ( | const DimensionedField< Type, volMesh > & | , |
const bool | allowUnmapped = false |
||
) | const |
Map volume internal (dimensioned) field Optional unmapped argument (currently unused)
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
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)
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)
|
static |
Name for exposed internal faces (default: oldInternalFaces)
Definition at line 164 of file fvMeshSubset.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.