Loading...
Searching...
No Matches
pointSmoother Class Referenceabstract

Abstract base class for point smoothing methods. Handles parallel communication via reset and average functions. More...

#include <pointSmoother.H>

Inheritance diagram for pointSmoother:

Public Member Functions

 TypeName ("pointSmoother")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, pointSmoother, dictionary,(const polyMesh &mesh, const dictionary &dict),(mesh, dict))
 pointSmoother (const polyMesh &mesh, const dictionary &)
 Construct from a dictionary and a point displacement field.
virtual ~pointSmoother ()=default
 Destructor.
const polyMeshmesh () const noexcept
 Access the mesh.
void update (const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes, pointVectorField &pointDisplacement, const bool correctBCs=true) const
 Update the point displacements and apply constraints.
void update (const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const polyMeshGeometry &meshGeometry, pointVectorField &pointDisplacement, const bool correctBCs=true) const
 Update the point displacements and apply constraints.
virtual void calculate (const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes, vectorField &pointDisplacement) const =0
 Calculate the point displacement.
virtual void calculate (const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const polyMeshGeometry &meshGeometry, vectorField &pointDisplacement) const
 Update the point displacements.
virtual tmp< scalarFieldfaceQuality (const pointField &points, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes) const
 Check element quality: 1 = best, 0 = invalid. (also negative?) Topology from mesh, point locations supplied. Move to motionSolver level?
virtual tmp< scalarFieldcellQuality (const pointField &points, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes) const
 Check element quality: 1 = best, 0 = invalid. Topology from mesh, point locations supplied. Move to motionSolver level?

Static Public Member Functions

static autoPtr< pointSmoother > New (const word &pointSmootherType, const polyMesh &mesh, const dictionary &dict)
 Construct given type.
static autoPtr< pointSmoother > New (const polyMesh &mesh, const dictionary &dict)
 Construct with type looked up from dictionary.

Protected Member Functions

bool isInternalOrProcessorFace (const label faceI) const
 Test if the given face is internal or on a processor boundary.
bitSet pointsToMove (const labelList &facesToMove, const bool moveInternalFaces) const
 Get list of the points to be moved.
template<class weightType>
void reset (const labelList &facesToMove, Field< weightType > &weights, vectorField &pointDisplacement, const bool resetInternalFaces=true) const
 Reset the relevant weights and displacements to zero.
template<class weightType>
void average (const labelList &facesToMove, Field< weightType > &weights, vectorField &pointDisplacement) const
 Average the displacements using the weights provided.

Detailed Description

Abstract base class for point smoothing methods. Handles parallel communication via reset and average functions.

Source files

Definition at line 50 of file pointSmoother.H.

Constructor & Destructor Documentation

◆ pointSmoother()

pointSmoother ( const polyMesh & mesh,
const dictionary & dict )

Construct from a dictionary and a point displacement field.

Definition at line 91 of file pointSmoother.C.

References dict, mesh(), and pp().

Here is the call graph for this function:

◆ ~pointSmoother()

virtual ~pointSmoother ( )
virtualdefault

Destructor.

Member Function Documentation

◆ isInternalOrProcessorFace()

bool isInternalOrProcessorFace ( const label faceI) const
protected

Test if the given face is internal or on a processor boundary.

Definition at line 36 of file pointSmoother.C.

References mesh(), and patchID.

Referenced by equipotentialPointSmoother::calculate(), geometricElementTransformPointSmoother::calculate(), laplacianPointSmoother::calculate(), and pointsToMove().

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

◆ pointsToMove()

Foam::bitSet pointsToMove ( const labelList & facesToMove,
const bool moveInternalFaces ) const
protected

Get list of the points to be moved.

Definition at line 61 of file pointSmoother.C.

References isInternalOrProcessorFace(), mesh, nPoints, bitSet::set(), syncTools::syncPointList(), and U.

Referenced by average(), laplacianConstraintPointSmoother::calculate(), and reset().

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

◆ reset()

template<class weightType>
void reset ( const labelList & facesToMove,
Field< weightType > & weights,
vectorField & pointDisplacement,
const bool resetInternalFaces = true ) const
protected

Reset the relevant weights and displacements to zero.

Definition at line 27 of file pointSmootherTemplates.C.

References pointsToMove(), and VectorSpace< Form, Cmpt, Ncmpts >::zero.

Here is the call graph for this function:

◆ average()

template<class weightType>
void average ( const labelList & facesToMove,
Field< weightType > & weights,
vectorField & pointDisplacement ) const
protected

Average the displacements using the weights provided.

Definition at line 49 of file pointSmootherTemplates.C.

References mesh, pointsToMove(), syncTools::syncPointList(), and VectorSpace< Form, Cmpt, Ncmpts >::zero.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "pointSmoother" )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr ,
pointSmoother ,
dictionary ,
(const polyMesh &mesh, const dictionary &dict) ,
(mesh, dict)  )

References dict, and mesh().

Here is the call graph for this function:

◆ New() [1/2]

Foam::autoPtr< Foam::pointSmoother > New ( const word & pointSmootherType,
const polyMesh & mesh,
const dictionary & dict )
static

Construct given type.

Definition at line 112 of file pointSmoother.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, Foam::Info, and mesh().

Here is the call graph for this function:

◆ New() [2/2]

Foam::autoPtr< Foam::pointSmoother > New ( const polyMesh & mesh,
const dictionary & dict )
static

Construct with type looked up from dictionary.

Definition at line 139 of file pointSmoother.C.

References dict, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, mesh(), and Foam::New().

Here is the call graph for this function:

◆ mesh()

◆ update() [1/2]

void update ( const labelList & facesToMove,
const pointField & oldPoints,
const pointField & currentPoints,
const pointField & faceCentres,
const vectorField & faceAreas,
const pointField & cellCentres,
const scalarField & cellVolumes,
pointVectorField & pointDisplacement,
const bool correctBCs = true ) const

Update the point displacements and apply constraints.

Definition at line 153 of file pointSmoother.C.

References calculate(), pointConstraints::constrainDisplacement(), mesh, MeshObject< pointMesh, UpdateableMeshObject, pointConstraints >::New(), and GeometricField< Type, PatchField, GeoMesh >::ref().

Here is the call graph for this function:

◆ update() [2/2]

void update ( const labelList & facesToMove,
const pointField & oldPoints,
const pointField & currentPoints,
const polyMeshGeometry & meshGeometry,
pointVectorField & pointDisplacement,
const bool correctBCs = true ) const
inline

Update the point displacements and apply constraints.

Definition at line 206 of file pointSmoother.H.

References polyMeshGeometry::cellCentres(), polyMeshGeometry::cellVolumes(), polyMeshGeometry::faceAreas(), polyMeshGeometry::faceCentres(), and update().

Here is the call graph for this function:

◆ calculate() [1/2]

virtual void calculate ( const labelList & facesToMove,
const pointField & oldPoints,
const pointField & currentPoints,
const pointField & faceCentres,
const vectorField & faceAreas,
const pointField & cellCentres,
const scalarField & cellVolumes,
vectorField & pointDisplacement ) const
pure virtual

Calculate the point displacement.

Implemented in equipotentialPointSmoother, geometricElementTransformPointSmoother, laplacianConstraintPointSmoother, and laplacianPointSmoother.

Referenced by calculate(), and update().

Here is the caller graph for this function:

◆ calculate() [2/2]

virtual void calculate ( const labelList & facesToMove,
const pointField & oldPoints,
const pointField & currentPoints,
const polyMeshGeometry & meshGeometry,
vectorField & pointDisplacement ) const
inlinevirtual

Update the point displacements.

Definition at line 248 of file pointSmoother.H.

References calculate(), polyMeshGeometry::cellCentres(), polyMeshGeometry::cellVolumes(), polyMeshGeometry::faceAreas(), and polyMeshGeometry::faceCentres().

Here is the call graph for this function:

◆ faceQuality()

Foam::tmp< Foam::scalarField > faceQuality ( const pointField & points,
const pointField & faceCentres,
const vectorField & faceAreas,
const pointField & cellCentres,
const scalarField & cellVolumes ) const
virtual

Check element quality: 1 = best, 0 = invalid. (also negative?) Topology from mesh, point locations supplied. Move to motionSolver level?

Definition at line 198 of file pointSmoother.C.

References polyMeshTools::faceOrthogonality(), mesh, and points.

Referenced by cellQuality().

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

◆ cellQuality()

Foam::tmp< Foam::scalarField > cellQuality ( const pointField & points,
const pointField & faceCentres,
const vectorField & faceAreas,
const pointField & cellCentres,
const scalarField & cellVolumes ) const
virtual

Check element quality: 1 = best, 0 = invalid. Topology from mesh, point locations supplied. Move to motionSolver level?

Definition at line 229 of file pointSmoother.C.

References polyMesh::faceNeighbour(), polyMesh::faceOwner(), faceQuality(), forAll, mesh, Foam::min(), tmp< T >::New(), primitiveMesh::nInternalFaces(), and points.

Here is the call graph for this function:

The documentation for this class was generated from the following files:
  • src/dynamicMesh/motionSolvers/displacement/pointSmoothing/pointSmoothers/pointSmoother/pointSmoother.H
  • src/dynamicMesh/motionSolvers/displacement/pointSmoothing/pointSmoothers/pointSmoother/pointSmoother.C
  • src/dynamicMesh/motionSolvers/displacement/pointSmoothing/pointSmoothers/pointSmoother/pointSmootherTemplates.C