cellCellStencil Class Referenceabstract

Calculation of interpolation stencils. More...

Inheritance diagram for cellCellStencil:
Collaboration diagram for cellCellStencil:

Public Types

enum  patchCellType { OTHER = 0, PATCH = 1, OVERSET = 2 }
 
enum  cellType {
  CALCULATED = 0, INTERPOLATED = 1, HOLE = 2, SPECIAL = 3,
  POROUS = 4
}
 

Public Member Functions

 TypeName ("cellCellStencil")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, cellCellStencil, mesh,(const fvMesh &mesh, const dictionary &dict, const bool update),(mesh, dict, update))
 
 cellCellStencil (const fvMesh &)
 Construct from fvMesh. More...
 
virtual ~cellCellStencil ()
 Destructor. More...
 
virtual bool update ()=0
 Update stencils. Return false if nothing changed. More...
 
virtual const labelUListcellTypes () const =0
 Return the cell type list. More...
 
virtual const labelUListinterpolationCells () const =0
 Indices of interpolated cells. More...
 
virtual const mapDistributecellInterpolationMap () const =0
 Return a communication schedule. More...
 
virtual const labelListListcellStencil () const =0
 Per interpolated cell the neighbour cells (in terms of slots as. More...
 
virtual const List< scalarList > & cellInterpolationWeights () const =0
 Weights for cellStencil. More...
 
virtual const scalarListcellInterpolationWeight () const =0
 Per interpolated cell the interpolation factor. (0 = use. More...
 
virtual void stencilWeights (const point &sample, const pointList &donorCcs, scalarList &weights) const =0
 Calculate weights for a single acceptor. More...
 
virtual const wordHashSetnonInterpolatedFields () const
 Return the names of any (stencil or mesh specific) fields that. More...
 
virtual wordHashSetnonInterpolatedFields ()
 Return non-const non-interpolating fields. More...
 
bool localStencil (const labelUList &) const
 Helper: is stencil fully local. More...
 
const labelIOListzoneID () const
 Helper: get reference to registered zoneID. Loads volScalarField. More...
 
template<class T >
void interpolate (const fvMesh &mesh, Field< T > &psi) const
 Explicit interpolation of acceptor cells from donor cells. Looks up cellCellStencil. More...
 
template<class GeoField >
void interpolate (GeoField &psi) const
 Explicit interpolation of acceptor cells from donor cells with boundary condition handling. More...
 
template<class GeoField >
void interpolate (const fvMesh &mesh, const wordHashSet &suppressed) const
 Explicit interpolation of all registered fields. No boundary conditions. Excludes selected fields (and their old-time fields) More...
 
void walkFront (const globalIndex &globalCells, const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight, const scalarList &compactCellVol, const labelListList &compactStencil, const labelList &zoneID, const label holeLayers, const label useLayer) const
 Surround holes with layer(s) of interpolated cells. More...
 
void setUpFront (const labelList &allCellTypes, bitSet &isFront) const
 Set up front using allCellTypes. More...
 
void setUpFrontOnOversetPatch (const labelList &allCellTypes, bitSet &isFront) const
 Set up front on overset patches. More...
 
void seedCell (const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
 Seed faces of cell with wantedFraction (if higher than current) More...
 
InfoProxy< cellCellStencilinfo () const noexcept
 Return info proxy, used to print stencil information to a stream. More...
 
template<class Type >
Foam::tmp< Foam::volScalarFieldcreateField (const fvMesh &mesh, const word &name, const UList< Type > &psi)
 

Static Public Member Functions

static autoPtr< cellCellStencilNew (const fvMesh &, const dictionary &dict, const bool update=true)
 New function which constructs and returns pointer to a. More...
 
static const labelIOListzoneID (const fvMesh &)
 Helper: get reference to registered zoneID. Loads volScalarField. More...
 
static labelList count (const label size, const labelUList &lst)
 Count occurrences (in parallel) More...
 
static void globalCellCells (const globalIndex &gi, const polyMesh &mesh, const boolList &isValidDonor, const labelList &selectedCells, labelListList &cellCells, pointListList &cellCellCentres)
 Helper: create cell-cell addressing in global numbering. More...
 
template<class T >
static void interpolate (Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
 Interpolation of acceptor cells from donor cells. More...
 
template<class GeoField , class SuppressBC >
static void correctBoundaryConditions (GeoField &psi)
 Version of correctBoundaryConditions that excludes 'overset' bcs. More...
 

Protected Member Functions

void suppressMotionFields ()
 Helper: populate nonInterpolatedFields_ with motion solver. More...
 

Static Protected Member Functions

template<class Type >
static tmp< volScalarFieldcreateField (const fvMesh &mesh, const word &name, const UList< Type > &)
 Helper: create volScalarField for postprocessing. More...
 
static word baseName (const word &name)
 Helper: strip off trailing _0. More...
 

Protected Attributes

const fvMeshmesh_
 Reference to the mesh. More...
 
wordHashSet nonInterpolatedFields_
 Set of fields that should not be interpolated. More...
 
const dictionary dict_
 Dictionary of motion control parameters. More...
 

Static Protected Attributes

static const Enum< cellTypecellTypeNames_
 Mode type names. More...
 

Detailed Description

Calculation of interpolation stencils.

Looks up zoneID labelIOList to give the zoning. Wrapped in MeshObject as cellCellStencilObject. Kept separate so meshes can implement more clever methods (e.g. solid body motion does not require full recalculation)

Source files

Definition at line 63 of file cellCellStencil.H.

Member Enumeration Documentation

◆ patchCellType

Enumerator
OTHER 
PATCH 
OVERSET 

Definition at line 67 of file cellCellStencil.H.

◆ cellType

enum cellType
Enumerator
CALCULATED 
INTERPOLATED 
HOLE 
SPECIAL 
POROUS 

Definition at line 74 of file cellCellStencil.H.

Constructor & Destructor Documentation

◆ cellCellStencil()

cellCellStencil ( const fvMesh mesh)

Construct from fvMesh.

Definition at line 52 of file cellCellStencil.C.

◆ ~cellCellStencil()

~cellCellStencil ( )
virtual

Destructor.

Definition at line 89 of file cellCellStencil.C.

Member Function Documentation

◆ createField() [1/2]

static tmp<volScalarField> createField ( const fvMesh mesh,
const word name,
const UList< Type > &   
)
staticprotected

Helper: create volScalarField for postprocessing.

◆ baseName()

Foam::word baseName ( const word name)
staticprotected

Helper: strip off trailing _0.

Definition at line 95 of file cellCellStencil.C.

References string::ends_with(), and Foam::name().

Here is the call graph for this function:

◆ suppressMotionFields()

void suppressMotionFields ( )
protected

Helper: populate nonInterpolatedFields_ with motion solver.

fields

Definition at line 106 of file cellCellStencil.C.

◆ TypeName()

TypeName ( "cellCellStencil"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
cellCellStencil  ,
mesh  ,
(const fvMesh &mesh, const dictionary &dict, const bool update ,
(mesh, dict, update  
)

◆ New()

Foam::autoPtr< Foam::cellCellStencil > New ( const fvMesh mesh,
const dictionary dict,
const bool  update = true 
)
static

New function which constructs and returns pointer to a.

cellCellStencil

Definition at line 60 of file cellCellStencil.C.

References DebugInFunction, dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::get(), mesh, and update().

Here is the call graph for this function:

◆ update()

virtual bool update ( )
pure virtual

Update stencils. Return false if nothing changed.

Implemented in inverseDistance, cellVolumeWeight, trackingInverseDistance, and cellCellStencilObject.

◆ cellTypes()

virtual const labelUList& cellTypes ( ) const
pure virtual

Return the cell type list.

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ interpolationCells()

virtual const labelUList& interpolationCells ( ) const
pure virtual

Indices of interpolated cells.

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ cellInterpolationMap()

virtual const mapDistribute& cellInterpolationMap ( ) const
pure virtual

Return a communication schedule.

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ cellStencil()

virtual const labelListList& cellStencil ( ) const
pure virtual

Per interpolated cell the neighbour cells (in terms of slots as.

constructed by above cellInterpolationMap) to interpolate

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ cellInterpolationWeights()

virtual const List<scalarList>& cellInterpolationWeights ( ) const
pure virtual

Weights for cellStencil.

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ cellInterpolationWeight()

virtual const scalarList& cellInterpolationWeight ( ) const
pure virtual

Per interpolated cell the interpolation factor. (0 = use.

calculated, 1 = use interpolated)

Implemented in inverseDistance, cellVolumeWeight, and cellCellStencilObject.

◆ stencilWeights()

virtual void stencilWeights ( const point sample,
const pointList donorCcs,
scalarList weights 
) const
pure virtual

Calculate weights for a single acceptor.

Implemented in inverseDistance, cellVolumeWeight, cellCellStencilObject, and leastSquares.

◆ nonInterpolatedFields() [1/2]

const Foam::wordHashSet & nonInterpolatedFields ( ) const
virtual

Return the names of any (stencil or mesh specific) fields that.

should not be interpolated

Reimplemented in cellCellStencilObject.

Definition at line 197 of file cellCellStencil.C.

◆ nonInterpolatedFields() [2/2]

Foam::wordHashSet & nonInterpolatedFields ( )
virtual

Return non-const non-interpolating fields.

Definition at line 203 of file cellCellStencil.C.

◆ localStencil()

bool localStencil ( const labelUList slots) const

Helper: is stencil fully local.

Definition at line 209 of file cellCellStencil.C.

References forAll.

◆ zoneID() [1/2]

const Foam::labelIOList & zoneID ( const fvMesh mesh)
static

Helper: get reference to registered zoneID. Loads volScalarField.

if not registered.

Definition at line 133 of file cellCellStencil.C.

References polyMesh::dbDir(), polyMesh::facesInstance(), Time::findInstance(), forAll, objectRegistry::getObjectPtr(), mesh, polyMesh::meshSubDir, IOobjectOption::MUST_READ, primitiveMesh::nCells(), IOobjectOption::NO_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, IOobjectOption::REGISTER, regIOobject::store(), and fvMesh::time().

Referenced by oversetFvMeshBase::writeObject().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ zoneID() [2/2]

const labelIOList& zoneID ( ) const
inline

Helper: get reference to registered zoneID. Loads volScalarField.

if not registered.

Definition at line 280 of file cellCellStencil.H.

References cellCellStencil::mesh_.

◆ count()

Foam::labelList count ( const label  size,
const labelUList lst 
)
static

Count occurrences (in parallel)

Definition at line 182 of file cellCellStencil.C.

References Foam::BitOps::count(), forAll, Pstream::listCombineGather(), and Foam::Zero.

Here is the call graph for this function:

◆ globalCellCells()

void globalCellCells ( const globalIndex gi,
const polyMesh mesh,
const boolList isValidDonor,
const labelList selectedCells,
labelListList cellCells,
pointListList cellCellCentres 
)
static

◆ interpolate() [1/4]

void interpolate ( Field< T > &  psi,
const fvMesh mesh,
const cellCellStencil overlap,
const List< scalarList > &  wghts 
)
static

Interpolation of acceptor cells from donor cells.

Definition at line 27 of file cellCellStencilTemplates.C.

References cellIDs, Foam::exit(), f(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, overlap, psi, s, UList< T >::size(), and T.

Referenced by cellCellStencil::interpolate(), and oversetFvMeshBase::interpolateFields().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interpolate() [2/4]

void interpolate ( const fvMesh mesh,
Field< T > &  psi 
) const

Explicit interpolation of acceptor cells from donor cells. Looks up cellCellStencil.

Definition at line 77 of file cellCellStencilTemplates.C.

References cellCellStencil::interpolate(), mesh, overlap, and psi.

Here is the call graph for this function:

◆ interpolate() [3/4]

void interpolate ( GeoField &  psi) const

Explicit interpolation of acceptor cells from donor cells with boundary condition handling.

Definition at line 92 of file cellCellStencilTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and psi.

Here is the call graph for this function:

◆ interpolate() [4/4]

void interpolate ( const fvMesh mesh,
const wordHashSet suppressed 
) const

Explicit interpolation of all registered fields. No boundary conditions. Excludes selected fields (and their old-time fields)

Definition at line 105 of file cellCellStencilTemplates.C.

References cellCellStencilObject::cellInterpolationWeights(), objectRegistry::csorted(), Foam::ensightOutput::debug, Foam::endl(), field(), fld, HashTable< T, Key, Hash >::found(), Foam::interpolate(), mesh, Foam::name(), overlap, Foam::Pout, and fvMesh::thisDb().

Here is the call graph for this function:

◆ correctBoundaryConditions()

void correctBoundaryConditions ( GeoField &  psi)
static

Version of correctBoundaryConditions that excludes 'overset' bcs.

Definition at line 181 of file cellCellStencilTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), forAll, UPstream::nonBlocking, UPstream::nRequests(), UPstream::parRun(), psi, and UPstream::waitRequests().

Referenced by oversetFvMeshBase::writeObject().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ walkFront()

void walkFront ( const globalIndex globalCells,
const scalar  layerRelax,
const labelListList allStencil,
labelList allCellTypes,
scalarField allWeight,
const scalarList compactCellVol,
const labelListList compactStencil,
const labelList zoneID,
const label  holeLayers,
const label  useLayer 
) const

◆ setUpFront()

void setUpFront ( const labelList allCellTypes,
bitSet isFront 
) const

Set up front using allCellTypes.

Definition at line 360 of file cellCellStencil.C.

References bitSet::set(), and syncTools::swapBoundaryCellList().

Here is the call graph for this function:

◆ setUpFrontOnOversetPatch()

void setUpFrontOnOversetPatch ( const labelList allCellTypes,
bitSet isFront 
) const

Set up front on overset patches.

Definition at line 408 of file cellCellStencil.C.

References fvBoundaryMesh::faceCells(), forAll, and bitSet::set().

Here is the call graph for this function:

◆ seedCell()

void seedCell ( const label  cellI,
const scalar  wantedFraction,
bitSet isFront,
scalarField fraction 
) const

Seed faces of cell with wantedFraction (if higher than current)

Definition at line 339 of file cellCellStencil.C.

References forAll, and bitSet::set().

Here is the call graph for this function:

◆ info()

InfoProxy<cellCellStencil> info ( ) const
inlinenoexcept

Return info proxy, used to print stencil information to a stream.

Definition at line 400 of file cellCellStencil.H.

◆ createField() [2/2]

Foam::tmp<Foam::volScalarField> createField ( const fvMesh mesh,
const word name,
const UList< Type > &  psi 
)

Member Data Documentation

◆ cellTypeNames_

const Foam::Enum< Foam::cellCellStencil::cellType > cellTypeNames_
staticprotected

Mode type names.

Definition at line 91 of file cellCellStencil.H.

◆ mesh_

const fvMesh& mesh_
protected

Reference to the mesh.

Definition at line 96 of file cellCellStencil.H.

Referenced by cellCellStencil::zoneID().

◆ nonInterpolatedFields_

wordHashSet nonInterpolatedFields_
protected

Set of fields that should not be interpolated.

Definition at line 101 of file cellCellStencil.H.

◆ dict_

const dictionary dict_
protected

Dictionary of motion control parameters.

Definition at line 106 of file cellCellStencil.H.


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