fvMeshSubsetProxy Class Reference

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 fvMeshbaseMesh () const noexcept
 The entire base mesh. More...
 
const fvMeshSubsetsubsetter () const noexcept
 The mesh subsetter. More...
 
fvMeshSubsetsubsetter () noexcept
 The mesh subsetter. More...
 
bool useSubMesh () const noexcept
 True if sub-mesh should be used. More...
 
const fvMeshmesh () const
 Access either base-mesh or sub-mesh. More...
 
const wordname () const noexcept
 The associated (set or zone) name if any. More...
 
const bitSetselectedCells () 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 > > 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. 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...
 

Detailed Description

Simple proxy for holding a mesh, or mesh-subset. The subMeshes are currently limited to cellSet or cellZone definitions.

Source files

Definition at line 55 of file fvMeshSubsetProxy.H.

Member Enumeration Documentation

◆ subsetType

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 62 of file fvMeshSubsetProxy.H.

Constructor & Destructor Documentation

◆ fvMeshSubsetProxy() [1/4]

fvMeshSubsetProxy ( fvMesh baseMesh)
explicit

Construct a pass-through proxy. No correct() invoked or required.

Definition at line 47 of file fvMeshSubsetProxy.C.

◆ fvMeshSubsetProxy() [2/4]

fvMeshSubsetProxy ( fvMesh baseMesh,
const subsetType  type,
const word selectionName,
label  exposedPatchId = -1 
)

Construct from components.

Definition at line 59 of file fvMeshSubsetProxy.C.

References correct(), UList< T >::first(), List< T >::resize(), fvMeshSubsetProxy::SET, fvMeshSubsetProxy::ZONE, and fvMeshSubsetProxy::ZONES.

Here is the call graph for this function:

◆ fvMeshSubsetProxy() [3/4]

fvMeshSubsetProxy ( fvMesh baseMesh,
const wordRes zoneNames,
label  exposedPatchId = -1 
)

Construct from components. The subsetType is ZONES.

Definition at line 90 of file fvMeshSubsetProxy.C.

References correct().

Here is the call graph for this function:

◆ fvMeshSubsetProxy() [4/4]

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

Here is the call graph for this function:

Member Function Documentation

◆ baseMesh()

const fvMesh & baseMesh ( ) const
inlinenoexcept

The entire base mesh.

Definition at line 147 of file fvMeshSubsetProxy.H.

Referenced by Foam::writeAllPointFields().

Here is the caller graph for this function:

◆ subsetter() [1/2]

const fvMeshSubset & subsetter ( ) const
inlinenoexcept

The mesh subsetter.

Definition at line 153 of file fvMeshSubsetProxy.H.

◆ subsetter() [2/2]

fvMeshSubset & subsetter ( )
inlinenoexcept

The mesh subsetter.

Definition at line 159 of file fvMeshSubsetProxy.H.

◆ useSubMesh()

bool useSubMesh ( ) const
inlinenoexcept

True if sub-mesh should be used.

Definition at line 165 of file fvMeshSubsetProxy.H.

References fvMeshSubsetProxy::NONE.

Referenced by fvMeshSubsetProxy::mesh(), and Foam::writePointField().

Here is the caller graph for this function:

◆ mesh()

const fvMesh & mesh ( ) const
inline

Access either base-mesh or sub-mesh.

Definition at line 171 of file fvMeshSubsetProxy.H.

References fvMeshSubset::subMesh(), and fvMeshSubsetProxy::useSubMesh().

Here is the call graph for this function:

◆ name()

const word & name ( ) const
inlinenoexcept

The associated (set or zone) name if any.

Definition at line 182 of file fvMeshSubsetProxy.H.

◆ selectedCells()

const bitSet & selectedCells ( ) const
inlinenoexcept

The current cell selection, when subsetting is active.

Definition at line 188 of file fvMeshSubsetProxy.H.

◆ resetZones()

void resetZones ( const wordRes zoneNames)

Define the zones selection, subset the mesh accordingly.

Definition at line 130 of file fvMeshSubsetProxy.C.

References correct(), and UList< T >::empty().

Here is the call graph for this function:

◆ correct()

bool correct ( bool  verbose = false)

Update of mesh subset.

Return true if the subset changed from previous call.

Definition at line 143 of file fvMeshSubsetProxy.C.

References Foam::endl(), Foam::flatOutput(), Foam::Info, PackedList< Width >::resize(), Foam::returnReduce(), bitSet::set(), and bitSet::transfer().

Here is the call graph for this function:

◆ readUpdate()

Read mesh. Correct on topo-change.

Definition at line 207 of file fvMeshSubsetProxy.C.

References correct(), polyMesh::POINTS_MOVED, polyMesh::TOPO_CHANGE, and polyMesh::TOPO_PATCH_CHANGE.

Here is the call graph for this function:

◆ zeroGradientField() [1/2]

static tmp< GeometricField< Type, fvPatchField, volMesh > > zeroGradientField ( const DimensionedField< Type, volMesh > &  df)
static

Construct volField (with zeroGradient) from an internal field.

◆ interpolateInternal() [1/8]

static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolateInternal ( const fvMeshSubset subsetter,
const DimensionedField< Type, volMesh > &  df 
)
static

Convert an internal field to a volume field (with zeroGradient)

◆ interpolateInternal() [2/8]

static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolateInternal ( const fvMeshSubset subsetter,
const tmp< DimensionedField< Type, volMesh > > &  tdf 
)
static

Convert an internal field to a volume field (with zeroGradient)

Currently no proper memory reuse

◆ interpolate() [1/8]

static tmp< GeoField > interpolate ( const fvMeshSubset subsetter,
const GeoField &  fld 
)
static

Wrapper for field or the subsetted field.

Pass through or forward to fvMeshSubset::interpolate()

Referenced by Foam::writePointField().

Here is the caller graph for this function:

◆ interpolate() [2/8]

static tmp< GeoField > interpolate ( const fvMeshSubset subsetter,
const tmp< GeoField > &  fld 
)
static

Wrapper for field or the subsetted field.

Pass through or forward to fvMeshSubset::interpolate()

◆ interpolateInternal() [3/8]

tmp< GeometricField< Type, fvPatchField, volMesh > > interpolateInternal ( const DimensionedField< Type, volMesh > &  df) const

Convert an internal field to a volume field (with zeroGradient)

◆ interpolateInternal() [4/8]

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

◆ interpolate() [3/8]

tmp< GeoField > interpolate ( const GeoField &  fld) const

Wrapper for field or the subsetted field.

Pass through or forward to fvMeshSubset::interpolate()

◆ interpolate() [4/8]

tmp< GeoField > interpolate ( const tmp< GeoField > &  fld) const

Wrapper for field or the subsetted field.

Pass through or forward to fvMeshSubset::interpolate()

◆ zeroGradientField() [2/2]

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > zeroGradientField ( const DimensionedField< Type, volMesh > &  df)

Definition at line 39 of file fvMeshSubsetProxyTemplates.C.

References DimensionedField< Type, GeoMesh >::dimensions(), io(), DimensionedField< Type, GeoMesh >::mesh(), Foam::New(), IOobject::NO_READ, IOobject::NO_WRITE, DimensionedField< Type, GeoMesh >::oriented(), IOobject::readOpt(), IOobject::registerObject(), IOobject::writeOpt(), and Foam::Zero.

Here is the call graph for this function:

◆ interpolateInternal() [5/8]

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal ( const fvMeshSubset subsetter,
const DimensionedField< Type, volMesh > &  df 
)

Definition at line 66 of file fvMeshSubsetProxyTemplates.C.

References fvMeshSubset::hasSubMesh(), and Foam::interpolate().

Here is the call graph for this function:

◆ interpolateInternal() [6/8]

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal ( const fvMeshSubset subsetter,
const tmp< DimensionedField< Type, volMesh > > &  tdf 
)

Definition at line 85 of file fvMeshSubsetProxyTemplates.C.

References fvMeshSubset::hasSubMesh(), and Foam::interpolate().

Here is the call graph for this function:

◆ interpolate() [5/8]

Foam::tmp< GeoField > interpolate ( const fvMeshSubset subsetter,
const GeoField &  fld 
)

Definition at line 120 of file fvMeshSubsetProxyTemplates.C.

References fld(), fvMeshSubset::hasSubMesh(), and fvMeshSubset::interpolate().

Here is the call graph for this function:

◆ interpolate() [6/8]

Foam::tmp< GeoField > interpolate ( const fvMeshSubset subsetter,
const tmp< GeoField > &  tfield 
)

Definition at line 141 of file fvMeshSubsetProxyTemplates.C.

References tmp< T >::clear(), fvMeshSubset::hasSubMesh(), Foam::interpolate(), and tmp< T >::valid().

Here is the call graph for this function:

◆ interpolateInternal() [7/8]

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal ( const DimensionedField< Type, volMesh > &  df) const

Definition at line 164 of file fvMeshSubsetProxyTemplates.C.

◆ interpolateInternal() [8/8]

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolateInternal ( const tmp< DimensionedField< Type, volMesh > > &  tdf) const

Definition at line 175 of file fvMeshSubsetProxyTemplates.C.

◆ interpolate() [7/8]

Foam::tmp< GeoField > interpolate ( const GeoField &  fld) const

Definition at line 186 of file fvMeshSubsetProxyTemplates.C.

References fld(), and Foam::interpolate().

Here is the call graph for this function:

◆ interpolate() [8/8]

Foam::tmp< GeoField > interpolate ( const tmp< GeoField > &  tfield) const

Definition at line 194 of file fvMeshSubsetProxyTemplates.C.

References Foam::interpolate().

Here is the call graph for this function:

The documentation for this class was generated from the following files: