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...
#include <fvMeshSubset.H>


Public Member Functions | |
| ClassName ("fvMeshSubset") | |
| fvMeshSubset (const fvMeshSubset &)=delete | |
| No copy construct. | |
| void | operator= (const fvMeshSubset &)=delete |
| No copy assignment. | |
| fvMeshSubset (const fvMesh &baseMesh) noexcept | |
| Construct using the entire mesh (no subset). | |
| fvMeshSubset (const fvMesh &baseMesh, Foam::zero) | |
| Construct a zero-sized subset mesh, non-processor patches only. | |
| fvMeshSubset (const fvMesh &baseMesh, const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true) | |
| Construct for a cell-subset of the given mesh. | |
| fvMeshSubset (const fvMesh &baseMesh, const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true) | |
| Construct for a cell-subset of the given mesh. | |
| fvMeshSubset (const fvMesh &baseMesh, const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true) | |
| Construct for a cell-subset of the given mesh. | |
| 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. | |
| const fvMesh & | baseMesh () const noexcept |
| Original mesh. | |
| const fvMesh & | mesh () const noexcept |
| Return baseMesh or subMesh, depending on the current state. | |
| bool | hasSubMesh () const noexcept |
| Have subMesh? | |
| const fvMesh & | subMesh () const |
| Return reference to subset mesh. | |
| fvMesh & | subMesh () |
| Return reference to subset mesh. | |
| const labelList & | pointMap () const |
| Return point map. | |
| const labelList & | faceMap () const |
| Return face map. | |
| const labelList & | faceFlipMap () const |
| Return face map with sign to encode flipped faces. | |
| const labelList & | cellMap () const |
| Return cell map. | |
| const labelList & | patchMap () const |
| Return patch map. | |
| const labelList & | pointPatchMap () const |
| Return point-patch map. Usually identical to patchMap except if additional patches are added to the pointMesh. | |
| void | clear () |
| Reset subMesh and all maps. | |
| void | reset () |
| Reset subMesh and all maps. Same as clear(). | |
| void | reset (Foam::zero) |
| Reset to a zero-sized subset mesh, non-processor patches only. | |
| void | reset (autoPtr< fvMesh > &&subMeshPtr, labelList &&pointMap, labelList &&faceMap, labelList &&cellMap, labelList &&patchMap) |
| Reset from components. | |
| void | reset (const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. | |
| void | reset (const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. | |
| void | reset (const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. | |
| 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. | |
| void | setCellSubset (const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. Same as reset(). | |
| void | setCellSubset (const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. Same as reset(). | |
| void | setCellSubset (const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true) |
| Use the specified subset of cells. Same as reset(). | |
| 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(). | |
| 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). | |
| template<class Type> | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const bool allowUnmapped=false) const |
| Map volume field. | |
| template<class Type> | |
| tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &, const bool allowUnmapped=false) const |
| Map surface field. | |
| template<class Type> | |
| tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &, const bool allowUnmapped=false) const |
| Map point field. | |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped) |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf, const bool allowUnmapped) const |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > | interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &vf, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap) |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > | interpolate (const GeometricField< Type, fvsPatchField, surfaceMesh > &sf, const bool allowUnmapped) const |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &vf, const pointMesh &sMesh, const labelUList &patchMap, const labelUList &pointMap) |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &sf, const bool allowUnmapped) const |
| template<class Type> | |
| Foam::tmp< Foam::DimensionedField< Type, Foam::volMesh > > | interpolate (const DimensionedField< Type, volMesh > &df, const fvMesh &sMesh, const labelUList &cellMap) |
| template<class Type> | |
| Foam::tmp< Foam::DimensionedField< Type, Foam::volMesh > > | interpolate (const DimensionedField< Type, volMesh > &df, const bool allowUnmapped) const |
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. | |
| 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. | |
| 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. | |
| 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. | |
Static Public Attributes | |
| static word | exposedPatchName |
| Name for exposed internal faces (default: oldInternalFaces). | |
Protected Member Functions | |
| bool | checkHasSubMesh () const |
| FatalError if subset has not been performed. | |
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.
|
delete |
No copy construct.
References fvMeshSubset().
Referenced by fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), and operator=().


|
explicitnoexcept |
Construct using the entire mesh (no subset).
Definition at line 400 of file fvMeshSubset.C.
References baseMesh(), and Foam::noexcept.

| fvMeshSubset | ( | const fvMesh & | baseMesh, |
| Foam::zero | ) |
Construct a zero-sized subset mesh, non-processor patches only.
Definition at line 406 of file fvMeshSubset.C.
References baseMesh(), fvMeshSubset(), and 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 414 of file fvMeshSubset.C.
References baseMesh(), fvMeshSubset(), patchID, and reset().

| 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 428 of file fvMeshSubset.C.
References baseMesh(), fvMeshSubset(), patchID, and reset().

| 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 442 of file fvMeshSubset.C.
References baseMesh(), fvMeshSubset(), patchID, and reset().

| 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 456 of file fvMeshSubset.C.
References baseMesh(), fvMeshSubset(), patchID, and reset().

|
protected |
FatalError if subset has not been performed.
Definition at line 112 of file fvMeshSubset.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and Foam::nl.
Referenced by cellMap(), faceMap(), patchMap(), pointMap(), pointPatchMap(), subMesh(), and subMesh().


| ClassName | ( | "fvMeshSubset" | ) |
|
delete |
No copy assignment.
References baseMesh(), cellMap(), clear(), faceFlipMap(), faceMap(), fvMeshSubset(), hasSubMesh(), mesh(), Foam::noexcept, patchID, patchMap(), pointMap(), pointPatchMap(), reset(), and subMesh().

|
inlinenoexcept |
Original mesh.
Definition at line 23 of file fvMeshSubsetI.H.
References Foam::noexcept.
Referenced by fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubsetter::getExposedFaces(), fvMeshSubsetter::getExposedFaces(), operator=(), reset(), reset(), reset(), reset(), fvMeshSubsetter::setCellSubset(), and fvMeshSubsetter::setCellSubset().

|
inlinenoexcept |
Return baseMesh or subMesh, depending on the current state.
Definition at line 29 of file fvMeshSubsetI.H.
References Foam::noexcept.
Referenced by operator=().

|
inlinenoexcept |
Have subMesh?
Definition at line 35 of file fvMeshSubsetI.H.
References Foam::noexcept.
Referenced by operator=().

|
inline |
Return reference to subset mesh.
Definition at line 41 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), marchingCells::initialise(), interpolate(), interpolate(), interpolate(), interpolate(), operator=(), structuredRenumber::renumber(), and momentumError::write().


|
inline |
Return reference to subset mesh.
Definition at line 49 of file fvMeshSubsetI.H.
References checkHasSubMesh().

|
inline |
Return point map.
Definition at line 57 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by fvMeshDistribute::distribute(), interpolate(), interpolate(), interpolate(), operator=(), and reset().


|
inline |
Return face map.
Definition at line 65 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by fvMeshDistribute::distribute(), marchingCells::initialise(), interpolate(), interpolate(), interpolate(), interpolate(), operator=(), and reset().


|
inline |
Return face map with sign to encode flipped faces.
Definition at line 73 of file fvMeshSubsetI.H.
Referenced by fvMeshDistribute::distribute(), and operator=().

|
inline |
Return cell map.
Definition at line 84 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by structuredDecomp::decompose(), fvMeshDistribute::distribute(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), operator=(), structuredRenumber::renumber(), and reset().


|
inline |
Return patch map.
Definition at line 92 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by fvMeshDistribute::distribute(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), interpolate(), operator=(), and reset().


|
inline |
Return point-patch map. Usually identical to patchMap except if additional patches are added to the pointMesh.
Definition at line 100 of file fvMeshSubsetI.H.
References checkHasSubMesh().
Referenced by interpolate(), and operator=().


| void clear | ( | ) |
Reset subMesh and all maps.
Definition at line 473 of file fvMeshSubset.C.
Referenced by operator=().

| void reset | ( | ) |
Reset subMesh and all maps. Same as clear().
Definition at line 485 of file fvMeshSubset.C.
References clear().
Referenced by fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), fvMeshSubset(), operator=(), setCellSubset(), setCellSubset(), setCellSubset(), and setCellSubset().


| void reset | ( | Foam::zero | ) |
Reset to a zero-sized subset mesh, non-processor patches only.
Definition at line 516 of file fvMeshSubset.C.
References clear(), DebugPout, Foam::endl(), forAll, Foam::identity(), Foam::isA(), MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), polyBoundaryMesh::nNonProcessor(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, List< label >::null(), PtrList< T >::set(), PtrList< T >::setSize(), and UPtrList< T >::size().

| 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 491 of file fvMeshSubset.C.
References cellMap(), clear(), faceMap(), patchMap(), and pointMap().

| void reset | ( | const bitSet & | selectedCells, |
| const label | patchID = -1, | ||
| const bool | syncPar = true ) |
Use the specified subset of cells.
Definition at line 625 of file fvMeshSubset.C.
References Foam::abort(), Pstream::allGatherList(), DynamicList< T, SizeMin >::append(), baseMesh(), UList< T >::begin(), clear(), DebugPout, Foam::endl(), exposedPatchName, Foam::FatalError, FatalErrorInFunction, polyBoundaryMesh::findPatchID(), forAll, Foam::identity(), Foam::inplaceReorder(), Foam::isA(), Foam::labelMax, Pstream::listReduce(), UPstream::myProcNo(), polyBoundaryMesh::names(), autoPtr< T >::New(), MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), newPointi, IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, UPstream::nProcs(), UPstream::parRun(), patchID, patchNames(), PtrList< T >::push_back(), Foam::reduce(), UPtrList< T >::reorder(), List< T >::resize(), PtrList< T >::resize(), PtrList< 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 1387 of file fvMeshSubset.C.
References baseMesh(), Foam::BitSetOps::create(), patchID, and reset().

| void reset | ( | const labelHashSet & | selectedCells, |
| const label | patchID = -1, | ||
| const bool | syncPar = true ) |
Use the specified subset of cells.
Definition at line 1403 of file fvMeshSubset.C.
References baseMesh(), Foam::BitSetOps::create(), patchID, and reset().

| 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 1419 of file fvMeshSubset.C.
References baseMesh(), Foam::BitSetOps::create(), patchID, and reset().

|
inline |
Use the specified subset of cells. Same as reset().
Definition at line 390 of file fvMeshSubset.H.
References patchID, reset(), and setCellSubset().
Referenced by setCellSubset().


|
inline |
Use the specified subset of cells. Same as reset().
Definition at line 403 of file fvMeshSubset.H.
References patchID, and reset().

|
inline |
Use the specified subset of cells. Same as reset().
Definition at line 416 of file fvMeshSubset.H.
References patchID, and reset().

|
inline |
Use the specified subset of cells. Same as reset().
Definition at line 429 of file fvMeshSubset.H.
References patchID, and reset().

|
static |
Map volume internal (dimensioned) field.
References cellMap().
Referenced by momentumError::calcMomentError(), and momentumError::divDevRhoReff().


|
static |
Map volume field.
Optionally allow unmapped faces not to produce a warning
References cellMap(), faceMap(), and patchMap().

|
static |
Map surface field.
Optionally negates value if flipping a face (from exposing an internal face)
References cellMap(), faceMap(), and patchMap().

|
static |
| 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)
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
| const fvMesh & | sMesh, | ||
| const labelUList & | patchMap, | ||
| const labelUList & | cellMap, | ||
| const labelUList & | faceMap, | ||
| const bool | allowUnmapped ) |
Definition at line 37 of file fvMeshSubsetTemplates.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), fvPatchFieldBase::calculatedType(), cellMap(), DimensionedField< Type, GeoMesh >::dimensions(), faceMap(), fld(), forAll, DirectFieldMapper< FieldMapperType >::hasUnmapped(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), fvPatchField< Type >::New(), Foam::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, DimensionedField< Type, volMesh >::null(), DimensionedField< Type, GeoMesh >::oriented(), patchMap(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), PtrList< T >::set(), fvPatch::size(), fvPatch::start(), fvMesh::time(), and Time::timeName().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
| const bool | allowUnmapped ) const |
Definition at line 194 of file fvMeshSubsetTemplates.C.
References cellMap(), Foam::faceMap(), Foam::interpolate(), patchMap(), and subMesh().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > interpolate | ( | const GeometricField< Type, fvsPatchField, surfaceMesh > & | vf, |
| const fvMesh & | sMesh, | ||
| const labelUList & | patchMap, | ||
| const labelUList & | cellMap, | ||
| const labelUList & | faceMap ) |
Definition at line 219 of file fvMeshSubsetTemplates.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), fvsPatchFieldBase::calculatedType(), cellMap(), DimensionedField< Type, GeoMesh >::dimensions(), faceMap(), forAll, GeometricField< Type, PatchField, GeoMesh >::internalField(), DimensionedField< Type, GeoMesh >::is_oriented(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), fvsPatchField< Type >::New(), Foam::New(), primitiveMesh::nInternalFaces(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, DimensionedField< Type, surfaceMesh >::null(), DimensionedField< Type, GeoMesh >::oriented(), fvPatch::patch(), patchMap(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), PtrList< T >::set(), fvPatch::size(), UList< T >::size(), fvPatch::start(), fvMesh::time(), Time::timeName(), and polyPatch::whichFace().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > interpolate | ( | const GeometricField< Type, fvsPatchField, surfaceMesh > & | sf, |
| const bool | allowUnmapped ) const |
Definition at line 383 of file fvMeshSubsetTemplates.C.
References cellMap(), Foam::faceMap(), Foam::interpolate(), patchMap(), and subMesh().

| Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | vf, |
| const pointMesh & | sMesh, | ||
| const labelUList & | patchMap, | ||
| const labelUList & | pointMap ) |
Definition at line 407 of file fvMeshSubsetTemplates.C.
References pointMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), pointPatchFieldBase::calculatedType(), HashTable< T, Key, Hash >::cfind(), DimensionedField< Type, GeoMesh >::dimensions(), forAll, Foam::invertToMap(), DimensionedField< Type, GeoMesh >::mesh(), pointPatch::meshPoints(), IOobject::name(), Foam::New(), pointPatchField< Type >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, DimensionedField< Type, pointMesh >::null(), DimensionedField< Type, GeoMesh >::oriented(), patchMap(), pointMap(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), PtrList< T >::set(), pointPatch::size(), pointMesh::thisDb(), pointMesh::time(), and Time::timeName().

| Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | sf, |
| const bool | allowUnmapped ) const |
Definition at line 537 of file fvMeshSubsetTemplates.C.
References Foam::interpolate(), MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), pointMap(), pointPatchMap(), and subMesh().

| Foam::tmp< Foam::DimensionedField< Type, Foam::volMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | df, |
| const fvMesh & | sMesh, | ||
| const labelUList & | cellMap ) |
Definition at line 560 of file fvMeshSubsetTemplates.C.
References cellMap(), DimensionedField< Type, GeoMesh >::dimensions(), IOobject::name(), Foam::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, DimensionedField< Type, GeoMesh >::oriented(), fvMesh::time(), and Time::timeName().

| Foam::tmp< Foam::DimensionedField< Type, Foam::volMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | df, |
| const bool | allowUnmapped ) const |
Definition at line 594 of file fvMeshSubsetTemplates.C.
References cellMap(), Foam::interpolate(), and subMesh().

|
static |
Name for exposed internal faces (default: oldInternalFaces).
Definition at line 162 of file fvMeshSubset.H.
Referenced by marchingCells::initialise(), and reset().