Abstract base class for point smoothing methods. Handles parallel communication via reset and average functions. More...
#include <pointSmoother.H>

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 polyMesh & | mesh () const noexcept |
| Access the mesh. | |
| void | update (const labelList &facesToMove, const pointField &oldPoints, const pointField ¤tPoints, 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 ¤tPoints, 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 ¤tPoints, 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 ¤tPoints, const polyMeshGeometry &meshGeometry, vectorField &pointDisplacement) const |
| Update the point displacements. | |
| virtual tmp< scalarField > | faceQuality (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< scalarField > | cellQuality (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. | |
Abstract base class for point smoothing methods. Handles parallel communication via reset and average functions.
Definition at line 50 of file pointSmoother.H.
| 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().

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


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


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

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

| TypeName | ( | "pointSmoother" | ) |
Runtime type information.
| declareRunTimeSelectionTable | ( | autoPtr | , |
| pointSmoother | , | ||
| dictionary | , | ||
| (const polyMesh &mesh, const dictionary &dict) | , | ||
| (mesh, 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().

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

|
inlinenoexcept |
Access the mesh.
Definition at line 182 of file pointSmoother.H.
References Foam::noexcept.
Referenced by geometricElementTransformPointSmoother::calculate(), geometricElementTransformPointSmoother::cellQuality(), declareRunTimeSelectionTable(), equipotentialPointSmoother::equipotentialPointSmoother(), geometricElementTransformPointSmoother::geometricElementTransformPointSmoother(), isInternalOrProcessorFace(), laplacianConstraintPointSmoother::laplacianConstraintPointSmoother(), laplacianPointSmoother::laplacianPointSmoother(), New(), New(), pointSmoother(), equipotentialPointSmoother::TypeName(), geometricElementTransformPointSmoother::TypeName(), laplacianConstraintPointSmoother::TypeName(), laplacianPointSmoother::TypeName(), and geometricElementTransformPointSmoother::~geometricElementTransformPointSmoother().

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

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

|
pure virtual |
Calculate the point displacement.
Implemented in equipotentialPointSmoother, geometricElementTransformPointSmoother, laplacianConstraintPointSmoother, and laplacianPointSmoother.
Referenced by calculate(), and update().

|
inlinevirtual |
Update the point displacements.
Definition at line 248 of file pointSmoother.H.
References calculate(), polyMeshGeometry::cellCentres(), polyMeshGeometry::cellVolumes(), polyMeshGeometry::faceAreas(), and polyMeshGeometry::faceCentres().

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


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