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 |
The entire base mesh. More... | |
const fvMeshSubset & | subsetter () const |
The mesh subsetter. More... | |
bool | useSubMesh () const |
Check if a sub-mesh is being used. More... | |
const fvMesh & | mesh () const |
Access either base-mesh or sub-mesh. More... | |
const word & | name () const |
The associated (set or zone) name if any. More... | |
const bitSet & | selectedCells () const |
The current cell selection, when subsetting is active. 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 55 of file fvMeshSubsetProxy.H.
enum subsetType |
Internal bookkeeping for subset type.
Enumerator | |
---|---|
NONE | |
SET | |
ZONE | |
ZONES |
Definition at line 62 of file fvMeshSubsetProxy.H.
|
explicit |
Construct a pass-through proxy. No correct() invoked or required.
Definition at line 35 of file fvMeshSubsetProxy.C.
References fvMeshSubsetProxy::correct(), and fvMeshSubsetProxy::useSubMesh().
fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
const subsetType | type, | ||
const word & | selectionName, | ||
label | exposedPatchId = -1 |
||
) |
Construct from components.
Definition at line 53 of file fvMeshSubsetProxy.C.
References correct().
fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
const wordRes & | zoneNames, | ||
label | exposedPatchId = -1 |
||
) |
Construct from components. The subsetType is ZONES.
Definition at line 87 of file fvMeshSubsetProxy.C.
References correct().
fvMeshSubsetProxy | ( | fvMesh & | baseMesh, |
wordRes && | zoneNames, | ||
label | exposedPatchId = -1 |
||
) |
Construct from components. The subsetType is ZONES.
Definition at line 109 of file fvMeshSubsetProxy.C.
References correct().
|
inline |
The entire base mesh.
Definition at line 144 of file fvMeshSubsetProxy.H.
Referenced by Foam::writeAllPointFields().
|
inline |
The mesh subsetter.
Definition at line 150 of file fvMeshSubsetProxy.H.
|
inline |
Check if a sub-mesh is being used.
Definition at line 156 of file fvMeshSubsetProxy.H.
References fvMeshSubsetProxy::NONE.
Referenced by fvMeshSubsetProxy::fvMeshSubsetProxy(), fvMeshSubsetProxy::mesh(), and Foam::writePointField().
|
inline |
Access either base-mesh or sub-mesh.
Definition at line 162 of file fvMeshSubsetProxy.H.
References fvMeshSubset::subMesh(), and fvMeshSubsetProxy::useSubMesh().
|
inline |
The associated (set or zone) name if any.
Definition at line 173 of file fvMeshSubsetProxy.H.
|
inline |
The current cell selection, when subsetting is active.
Definition at line 179 of file fvMeshSubsetProxy.H.
Update of mesh subset.
Return true if the subset changed from previous call.
Definition at line 132 of file fvMeshSubsetProxy.C.
References Foam::endl(), Foam::flatOutput(), Foam::Info, PackedList< Width >::resize(), Foam::returnReduce(), bitSet::set(), and bitSet::transfer().
Referenced by fvMeshSubsetProxy::fvMeshSubsetProxy().
Foam::polyMesh::readUpdateState readUpdate | ( | ) |
Read mesh. Correct on topo-change.
Definition at line 196 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 67 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 86 of file fvMeshSubsetProxyTemplates.C.
References fvMeshSubset::hasSubMesh(), and Foam::interpolate().
Foam::tmp<GeoField> interpolate | ( | const fvMeshSubset & | subsetter, |
const GeoField & | fld | ||
) |
Definition at line 121 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 142 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 165 of file fvMeshSubsetProxyTemplates.C.
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> > interpolateInternal | ( | const tmp< DimensionedField< Type, volMesh >> & | tdf | ) | const |
Definition at line 176 of file fvMeshSubsetProxyTemplates.C.
Foam::tmp<GeoField> interpolate | ( | const GeoField & | fld | ) | const |
Definition at line 186 of file fvMeshSubsetProxyTemplates.C.
References fld, and Foam::interpolate().
Definition at line 194 of file fvMeshSubsetProxyTemplates.C.
References Foam::interpolate().