Simple proxy for holding a mesh, or mesh-subset. The subMeshes are currently limited to cellSet or cellZone definitions. More...
#include <fvMeshSubsetProxy.H>
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. | |
| fvMeshSubsetProxy (fvMesh &baseMesh, const subsetType type, const word &selectionName, label exposedPatchId=-1) | |
| Construct from components. | |
| fvMeshSubsetProxy (fvMesh &baseMesh, const wordRes &zoneNames, label exposedPatchId=-1) | |
| Construct from components. The subsetType is ZONES. | |
| fvMeshSubsetProxy (fvMesh &baseMesh, wordRes &&zoneNames, label exposedPatchId=-1) | |
| Construct from components. The subsetType is ZONES. | |
| const fvMesh & | baseMesh () const noexcept |
| The entire base mesh. | |
| const fvMeshSubset & | subsetter () const noexcept |
| The mesh subsetter. | |
| fvMeshSubset & | subsetter () noexcept |
| The mesh subsetter. | |
| bool | useSubMesh () const noexcept |
| True if sub-mesh should be used. | |
| const fvMesh & | mesh () const |
| Access either base-mesh or sub-mesh. | |
| const word & | name () const noexcept |
| The associated (set or zone) name if any. | |
| const bitSet & | selectedCells () const noexcept |
| The current cell selection, when subsetting is active. | |
| void | resetZones (const wordRes &zoneNames) |
| Define the zones selection, subset the mesh accordingly. | |
| bool | correct (bool verbose=false) |
| Update of mesh subset. | |
| polyMesh::readUpdateState | readUpdate () |
| Read mesh. Correct on topo-change. | |
| 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). | |
| 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). | |
| template<class GeoField> | |
| tmp< GeoField > | interpolate (const GeoField &fld) const |
| Wrapper for field or the subsetted field. | |
| template<class GeoField> | |
| tmp< GeoField > | interpolate (const tmp< GeoField > &fld) const |
| Wrapper for field or the subsetted field. | |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | zeroGradientField (const DimensionedField< Type, volMesh > &df) |
| 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. | |
| 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). | |
| 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). | |
| template<class GeoField> | |
| static tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const GeoField &fld) |
| Wrapper for field or the subsetted field. | |
| template<class GeoField> | |
| static tmp< GeoField > | interpolate (const fvMeshSubset &subsetter, const tmp< GeoField > &fld) |
| Wrapper for field or the subsetted field. | |
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.
References baseMesh(), and NONE.

| fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
| const subsetType | type, | ||
| const word & | selectionName, | ||
| label | exposedPatchId = -1 ) |
Construct from components.
Definition at line 52 of file fvMeshSubsetProxy.C.
References baseMesh(), correct(), NONE, SET, ZONE, and ZONES.

Construct from components. The subsetType is ZONES.
Definition at line 83 of file fvMeshSubsetProxy.C.
References baseMesh(), correct(), and ZONES.

Construct from components. The subsetType is ZONES.
Definition at line 102 of file fvMeshSubsetProxy.C.
References baseMesh(), correct(), and ZONES.

|
inlinenoexcept |
The entire base mesh.
Definition at line 175 of file fvMeshSubsetProxy.H.
References Foam::noexcept.
Referenced by fvMeshSubsetProxy(), fvMeshSubsetProxy(), fvMeshSubsetProxy(), fvMeshSubsetProxy(), and Foam::writeAllPointFields().

|
inlinenoexcept |
The mesh subsetter.
Definition at line 183 of file fvMeshSubsetProxy.H.
References Foam::noexcept.
Referenced by interpolate(), interpolate(), interpolate(), interpolateInternal(), and interpolateInternal().

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

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

|
inlinenoexcept |
The associated (set or zone) name if any.
Definition at line 220 of file fvMeshSubsetProxy.H.
References Foam::noexcept.
|
inlinenoexcept |
The current cell selection, when subsetting is active.
Definition at line 228 of file fvMeshSubsetProxy.H.
References Foam::noexcept.
Referenced by correct().

| void resetZones | ( | const wordRes & | zoneNames | ) |
Define the zones selection, subset the mesh accordingly.
Definition at line 123 of file fvMeshSubsetProxy.C.
References correct(), UList< T >::empty(), and ZONES.

| 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, IOobjectOption::NO_REGISTER, NONE, Foam::returnReduceOr(), selectedCells(), SET, ZONE, and ZONES.

| 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.
Referenced by interpolateInternal(), and interpolateInternal().

|
static |
Convert an internal field to a volume field (with zeroGradient).
References subsetter().
Referenced by interpolateInternal(), and interpolateInternal().


|
static |
Convert an internal field to a volume field (with zeroGradient).
Currently no proper memory reuse
References subsetter().

|
static |
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
References fld(), and subsetter().
Referenced by Foam::writePointField().


|
static |
Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
References fld(), and subsetter().

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

Wrapper for field or the subsetted field.
Pass through or forward to fvMeshSubset::interpolate()
References fld().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > zeroGradientField | ( | const DimensionedField< Type, volMesh > & | df | ) |
Definition at line 32 of file fvMeshSubsetProxyTemplates.C.
References DimensionedField< Type, GeoMesh >::dimensions(), io, DimensionedField< Type, GeoMesh >::mesh(), Foam::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, DimensionedField< Type, GeoMesh >::oriented(), Foam::Zero, and fvPatchFieldBase::zeroGradientType().

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

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

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

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

| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal | ( | const DimensionedField< Type, volMesh > & | df | ) | const |
Definition at line 157 of file fvMeshSubsetProxyTemplates.C.
References interpolateInternal().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal | ( | const tmp< DimensionedField< Type, volMesh > > & | tdf | ) | const |
Definition at line 168 of file fvMeshSubsetProxyTemplates.C.
References interpolateInternal().

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