Simple proxy for holding a mesh, or mesh-subset. The subMeshes are currently limited to cellSet or cellZone definitions. More...
Public Types | |
| enum | subsetType { NONE, SET, ZONE, ZONES } |
| Internal bookkeeping for subset type. More... | |
Public Member Functions | |
| fvMeshSubsetProxy (fvMesh &baseMesh) | |
| Construct a pass-through proxy. No correct() invoked or required. More... | |
| fvMeshSubsetProxy (fvMesh &baseMesh, const subsetType type, const word &selectionName, label exposedPatchId=-1) | |
| Construct from components. More... | |
| fvMeshSubsetProxy (fvMesh &baseMesh, const wordRes &zoneNames, label exposedPatchId=-1) | |
| Construct from components. The subsetType is ZONES. More... | |
| fvMeshSubsetProxy (fvMesh &baseMesh, wordRes &&zoneNames, label exposedPatchId=-1) | |
| Construct from components. The subsetType is ZONES. More... | |
| const fvMesh & | baseMesh () const noexcept |
| The entire base mesh. More... | |
| const fvMeshSubset & | subsetter () const noexcept |
| The mesh subsetter. More... | |
| fvMeshSubset & | subsetter () noexcept |
| The mesh subsetter. More... | |
| bool | useSubMesh () const noexcept |
| True if sub-mesh should be used. More... | |
| const fvMesh & | mesh () const |
| Access either base-mesh or sub-mesh. More... | |
| const word & | name () const noexcept |
| The associated (set or zone) name if any. More... | |
| const bitSet & | selectedCells () const noexcept |
| The current cell selection, when subsetting is active. More... | |
| void | resetZones (const wordRes &zoneNames) |
| Define the zones selection, subset the mesh accordingly. More... | |
| bool | correct (bool verbose=false) |
| Update of mesh subset. More... | |
| polyMesh::readUpdateState | readUpdate () |
| Read mesh. Correct on topo-change. More... | |
| template<class Type > | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolateInternal (const DimensionedField< Type, volMesh > &df) const |
| Convert an internal field to a volume field (with zeroGradient) More... | |
| template<class Type > | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolateInternal (const tmp< DimensionedField< Type, volMesh >> &tdf) const |
| Convert an internal field to a volume field (with zeroGradient) More... | |
| template<class GeoField > | |
| tmp< GeoField > | interpolate (const GeoField &fld) const |
| Wrapper for field or the subsetted field. More... | |
| template<class GeoField > | |
| tmp< GeoField > | interpolate (const tmp< GeoField > &fld) const |
| Wrapper for field or the subsetted field. More... | |
| template<class Type > | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolateInternal (const fvMeshSubset &subsetter, const DimensionedField< Type, volMesh > &df) |
| template<class Type > | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolateInternal (const fvMeshSubset &subsetter, const tmp< DimensionedField< Type, volMesh >> &tdf) |
| template<class GeoField > | |
| Foam::tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const GeoField &fld) |
| template<class GeoField > | |
| Foam::tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const tmp< GeoField > &tfield) |
| template<class Type > | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolateInternal (const DimensionedField< Type, volMesh > &df) const |
| template<class Type > | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolateInternal (const tmp< DimensionedField< Type, volMesh >> &tdf) const |
| template<class GeoField > | |
| Foam::tmp< GeoField > | interpolate (const GeoField &fld) const |
| template<class GeoField > | |
| Foam::tmp< GeoField > | interpolate (const tmp< GeoField > &tfield) const |
Static Public Member Functions | |
| template<class Type > | |
| static tmp< GeometricField< Type, fvPatchField, volMesh > > | zeroGradientField (const DimensionedField< Type, volMesh > &df) |
| Construct volField (with zeroGradient) from an internal field. More... | |
| template<class Type > | |
| static tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolateInternal (const fvMeshSubset &subsetter, const DimensionedField< Type, volMesh > &df) |
| Convert an internal field to a volume field (with zeroGradient) More... | |
| template<class Type > | |
| static tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolateInternal (const fvMeshSubset &subsetter, const tmp< DimensionedField< Type, volMesh >> &tdf) |
| Convert an internal field to a volume field (with zeroGradient) More... | |
| template<class GeoField > | |
| static tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const GeoField &fld) |
| Wrapper for field or the subsetted field. More... | |
| template<class GeoField > | |
| static tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const tmp< GeoField > &fld) |
| Wrapper for field or the subsetted field. More... | |
Simple proxy for holding a mesh, or mesh-subset. The subMeshes are currently limited to cellSet or cellZone definitions.
Definition at line 51 of file fvMeshSubsetProxy.H.
| enum subsetType |
Internal bookkeeping for subset type.
| Enumerator | |
|---|---|
| NONE | No subset. |
| SET | Subset with a cellSet. |
| ZONE | Subset with a cellZone. |
| ZONES | Subset with multiple cellZones. |
Definition at line 60 of file fvMeshSubsetProxy.H.
|
explicit |
Construct a pass-through proxy. No correct() invoked or required.
Definition at line 40 of file fvMeshSubsetProxy.C.
| fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
| const subsetType | type, | ||
| const word & | selectionName, | ||
| label | exposedPatchId = -1 |
||
| ) |
| fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
| const wordRes & | zoneNames, | ||
| label | exposedPatchId = -1 |
||
| ) |
Construct from components. The subsetType is ZONES.
Definition at line 84 of file fvMeshSubsetProxy.C.
References correct.
| fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
| wordRes && | zoneNames, | ||
| label | exposedPatchId = -1 |
||
| ) |
Construct from components. The subsetType is ZONES.
Definition at line 103 of file fvMeshSubsetProxy.C.
References correct.
|
inlinenoexcept |
The entire base mesh.
Definition at line 175 of file fvMeshSubsetProxy.H.
Referenced by Foam::writeAllPointFields().

|
inlinenoexcept |
The mesh subsetter.
Definition at line 183 of file fvMeshSubsetProxy.H.
|
inlinenoexcept |
The mesh subsetter.
Definition at line 191 of file fvMeshSubsetProxy.H.
|
inlinenoexcept |
True if sub-mesh should be used.
Definition at line 199 of file fvMeshSubsetProxy.H.
Referenced by fvMeshSubsetProxy::mesh(), and Foam::writePointField().

|
inline |
Access either base-mesh or sub-mesh.
Definition at line 207 of file fvMeshSubsetProxy.H.
References fvMeshSubset::subMesh(), and fvMeshSubsetProxy::useSubMesh().

|
inlinenoexcept |
The associated (set or zone) name if any.
Definition at line 220 of file fvMeshSubsetProxy.H.
|
inlinenoexcept |
The current cell selection, when subsetting is active.
Definition at line 228 of file fvMeshSubsetProxy.H.
| void resetZones | ( | const wordRes & | zoneNames | ) |
Define the zones selection, subset the mesh accordingly.
Definition at line 123 of file fvMeshSubsetProxy.C.
References correct, and UList< T >::empty().

| bool correct | ( | bool | verbose = false | ) |
Update of mesh subset.
Return true if the subset changed from previous call.
Definition at line 136 of file fvMeshSubsetProxy.C.
References Foam::endl(), Foam::flatOutput(), Foam::Info, PackedList< Width >::resize(), Foam::returnReduceOr(), bitSet::set(), and bitSet::transfer().

| Foam::polyMesh::readUpdateState readUpdate | ( | ) |
Read mesh. Correct on topo-change.
Definition at line 200 of file fvMeshSubsetProxy.C.
References correct, polyMesh::POINTS_MOVED, polyMesh::TOPO_CHANGE, and polyMesh::TOPO_PATCH_CHANGE.
|
static |
Construct volField (with zeroGradient) from an internal field.
|
static |
Convert an internal field to a volume field (with zeroGradient)
|
static |
Convert an internal field to a volume field (with zeroGradient)
Currently no proper memory reuse
|
static |
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
Referenced by Foam::writePointField().

|
static |
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
| tmp<GeometricField<Type, fvPatchField, volMesh> > interpolateInternal | ( | const DimensionedField< Type, volMesh > & | df | ) | const |
Convert an internal field to a volume field (with zeroGradient)
| tmp<GeometricField<Type, fvPatchField, volMesh> > interpolateInternal | ( | const tmp< DimensionedField< Type, volMesh >> & | tdf | ) | const |
Convert an internal field to a volume field (with zeroGradient)
Currently no proper memory reuse
| tmp<GeoField> interpolate | ( | const GeoField & | fld | ) | const |
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
| Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> > interpolateInternal | ( | const fvMeshSubset & | subsetter, |
| const DimensionedField< Type, volMesh > & | df | ||
| ) |
Definition at line 60 of file fvMeshSubsetProxyTemplates.C.
References fvMeshSubset::hasSubMesh(), and Foam::interpolate().

| Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> > interpolateInternal | ( | const fvMeshSubset & | subsetter, |
| const tmp< DimensionedField< Type, volMesh >> & | tdf | ||
| ) |
Definition at line 79 of file fvMeshSubsetProxyTemplates.C.
References fvMeshSubset::hasSubMesh(), and Foam::interpolate().

| Foam::tmp<GeoField> interpolate | ( | const fvMeshSubset & | subsetter, |
| const GeoField & | fld | ||
| ) |
Definition at line 114 of file fvMeshSubsetProxyTemplates.C.
References fld, fvMeshSubset::hasSubMesh(), and fvMeshSubset::interpolate().

| Foam::tmp<GeoField> interpolate | ( | const fvMeshSubset & | subsetter, |
| const tmp< GeoField > & | tfield | ||
| ) |
Definition at line 135 of file fvMeshSubsetProxyTemplates.C.
References tmp< T >::clear(), fvMeshSubset::hasSubMesh(), Foam::interpolate(), and tmp< T >::valid().

| Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> > interpolateInternal | ( | const DimensionedField< Type, volMesh > & | df | ) | const |
Definition at line 158 of file fvMeshSubsetProxyTemplates.C.
| Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> > interpolateInternal | ( | const tmp< DimensionedField< Type, volMesh >> & | tdf | ) | const |
Definition at line 169 of file fvMeshSubsetProxyTemplates.C.
| Foam::tmp<GeoField> interpolate | ( | const GeoField & | fld | ) | const |
Definition at line 179 of file fvMeshSubsetProxyTemplates.C.
References fld, and Foam::interpolate().

Definition at line 187 of file fvMeshSubsetProxyTemplates.C.
References Foam::interpolate().

Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.