rsvs3D  0.0.0
Codes for the c++ implementation of the 3D RSVS
meshprocessing.hpp File Reference

Tools for the mathematical processing of meshes. More...

#include <vector>
#include "arraystructures.hpp"
#include "mesh.hpp"
#include "snake.hpp"

Go to the source code of this file.

Functions

void FlattenBoundaryFaces (mesh &meshin)
 
void TriangulateAllFaces (mesh &meshin)
 
std::vector< int > FindHolesInSnake (const snake &snakein, const HashedVector< int, int > &uncertainVert)
 
void PrepareSnakeForCFD (const snake &snakein, double distanceTol, mesh &meshgeom, std::vector< double > &holeCoords)
 Prepares the snake to be used for CFD, removes duplicate points and triangulates it. More...
 
HashedVector< int, int > GroupCloseSnaxels (const snake &snakein, double distTol)
 
void TestVertClose (int vertIndIn, std::vector< bool > &isSnaxDone, const mesh &meshin, double distTol, std::vector< int > &sameEdges)
 
HashedVector< int, int > GroupCloseVertices (const mesh &meshin, double distTol)
 
int FindVertexHole (int vertInd, const mesh &meshin, const std::vector< bool > &vertIn, const HashedVector< int, int > &uncertainVert, std::vector< bool > &vertExplored)
 
double DomInter (double x, double y1, double y2)
 
mesh BuildDomain (const std::array< double, 3 > &lowerB, const std::array< double, 3 > &upperB, double tolInner=0.0)
 
mesh BuildCutCellDomain (const std::array< double, 3 > &outerLowerB, const std::array< double, 3 > &outerUpperB, const std::array< double, 3 > &innerLowerB, const std::array< double, 3 > &innerUpperB, int nSteps, std::vector< int > &vertPerSubDomain)
 Builds a series of domains with different edge properties controlling the interpolation of the metric. More...
 
double PseudoSurfaceAngle (const mesh &meshin, const std::array< int, 2 > &surfInds)
 Calculates the pseudo surface angle. More...
 
std::vector< double > CalculateEdgeCurvature (const mesh &meshin)
 Calculates the angles between the surfaces connected at an edge. More...
 
std::vector< double > CalculateVertexCurvature (const mesh &meshin, int smoothingSteps)
 Calculates the vertex curvature. More...
 
std::vector< double > CalculateVertexMinEdgeLength (const mesh &meshin)
 Calculates the vertex minimum edge length. More...
 
std::vector< double > CalculateVertexMaxEdgeLength (const mesh &meshin)
 Calculates the vertex maximum edge length. More...
 
std::vector< double > CalculateVertexMeanEdgeLength (const mesh &meshin)
 Calculates the vertex mean edge length. More...
 
std::vector< double > CalculateEdgeLengths (const mesh &meshin)
 Calculates the edge lengths. More...
 
std::vector< double > CoordInVolume (const mesh &meshin)
 Generate a vector of coordinates of points probably inside volumes. More...
 
std::vector< double > VolumeCentroids (const mesh &meshin)
 Generate a vector of coordinates of points at the volume pseudo centroid. More...
 
std::vector< double > VolumeInternalLayers (const mesh &meshin, int nLayers)
 Returns points on edges between volume pseudo centroid and vertices. More...
 
std::vector< double > SurfaceCentroids (const mesh &meshin)
 Generate a vector of coordinates of points at the surfaces pseudo centroid. More...
 
std::vector< double > SurfaceInternalLayers (const mesh &meshin, int nLayers)
 Returns points on edges between surface pseudo centroid and vertices. More...
 
double CotanLaplacianWeight (const std::vector< double > &centre, const std::vector< double > &p_im1, const std::vector< double > &p_i, const std::vector< double > &p_ip1, coordvec &temp1, coordvec &temp2)
 
int VertexLaplacianVector (const mesh &meshin, const vert *vertsmooth, coordvec &lapVec, bool isOrdered=false)
 
coordvec VertexLaplacianVector (const mesh &meshin, int vertIndex)
 
std::array< double, 2 > IntersectLineSphere (const coordvec &lineVec, const coordvec &offset, double sphereRadius)
 Calculates the parametric position along a line at which it intersects a sphere. More...
 
std::array< double, 2 > IntersectLineSphere (const coordvec &lineVec, const std::vector< double > &linePoint, const coordvec &sphereCentre, double sphereRadius)
 
std::vector< double > MeshUnitNormals (const mesh &meshin)
 
std::vector< double > MeshLaplacians (const mesh &meshin)
 

Detailed Description

Tools for the mathematical processing of meshes.

Definition in file meshprocessing.hpp.

Function Documentation

◆ BuildCutCellDomain()

mesh BuildCutCellDomain ( const std::array< double, 3 > &  outerLowerB,
const std::array< double, 3 > &  outerUpperB,
const std::array< double, 3 > &  innerLowerB,
const std::array< double, 3 > &  innerUpperB,
int  nSteps,
std::vector< int > &  vertPerSubDomain 
)

Builds a series of domains with different edge properties controlling the interpolation of the metric.

Parameters
[in]outerLowerBThe outer lower b
[in]outerUpperBThe outer upper b
[in]innerLowerBThe inner lower b
[in]innerUpperBThe inner upper b
[in]nStepsThe steps
vertPerSubDomainThe vertical per sub domain
Returns
The cut cell domain.
 These also serve to avoid having a badly conditioned initial
triangulation with very small edges.
   `nSteps` is the number of total domains.
   0 will return an empty mesh, 1 will return a mesh of the inner bound
   2 will return inner and outer bounds,

                   {
                           {DomInter(x, innerLowerB[0], outerLowerB[0]),
                           DomInter(x, innerUpperB[0], outerUpperB[0])},
                           {DomInter(x, innerLowerB[1], outerLowerB[1]),
                           DomInter(x, innerUpperB[1], outerUpperB[1])},
                           {DomInter(x, innerLowerB[2], outerLowerB[2]),
                           DomInter(x, innerUpperB[2], outerUpperB[2])}
                   }

           meshtemp = BuildDomain(
                   {DomInter(x, innerLowerB[0], outerLowerB[0]),
                           DomInter(x, innerLowerB[1], outerLowerB[1]),
                           DomInter(x, innerLowerB[2], outerLowerB[2])},
                   {DomInter(x, innerUpperB[0], outerUpperB[0]),
                           DomInter(x, innerUpperB[1], outerUpperB[1]),
                           DomInter(x, innerUpperB[2], outerUpperB[2])}
                   );

Definition at line 512 of file meshprocessing.cpp.

◆ CalculateEdgeCurvature()

std::vector<double> CalculateEdgeCurvature ( const mesh meshin)

Calculates the angles between the surfaces connected at an edge.

To work the faces need have a common orientation

Parameters
[in]meshinThe input mesh
Returns
The edge angles.

Definition at line 612 of file meshprocessing.cpp.

◆ CalculateEdgeLengths()

std::vector<double> CalculateEdgeLengths ( const mesh meshin)

Calculates the edge lengths.

Parameters
[in]meshinThe meshin
Returns
The edge lengths.

Definition at line 773 of file meshprocessing.cpp.

◆ CalculateVertexCurvature()

std::vector<double> CalculateVertexCurvature ( const mesh meshin,
int  smoothingSteps 
)

Calculates the vertex curvature.

Parameters
[in]meshinThe input mesh
[in]smoothingStepsThe number of metric smoothing steps
Returns
The vertex curvature.

Definition at line 641 of file meshprocessing.cpp.

◆ CalculateVertexMaxEdgeLength()

std::vector<double> CalculateVertexMaxEdgeLength ( const mesh meshin)

Calculates the vertex maximum edge length.

Parameters
[in]meshinThe meshin
Returns
The vertex minimum edge length.

Definition at line 718 of file meshprocessing.cpp.

◆ CalculateVertexMeanEdgeLength()

std::vector<double> CalculateVertexMeanEdgeLength ( const mesh meshin)

Calculates the vertex mean edge length.

Parameters
[in]meshinThe meshin
Returns
The vertex mean edge length.

Definition at line 745 of file meshprocessing.cpp.

◆ CalculateVertexMinEdgeLength()

std::vector<double> CalculateVertexMinEdgeLength ( const mesh meshin)

Calculates the vertex minimum edge length.

Parameters
[in]meshinThe meshin
Returns
The vertex minimum edge length.

Definition at line 691 of file meshprocessing.cpp.

◆ CoordInVolume()

std::vector<double> CoordInVolume ( const mesh meshin)

Generate a vector of coordinates of points probably inside volumes.

Parameters
[in]meshinThe input mesh
Returns
vector of coordinates with one coordinate inside each volume

Definition at line 794 of file meshprocessing.cpp.

◆ IntersectLineSphere()

std::array<double, 2> IntersectLineSphere ( const coordvec lineVec,
const coordvec offset,
double  sphereRadius 
)

Calculates the parametric position along a line at which it intersects a sphere.

The return value of this function when multiplied by the line direction and added to the reference point gives the coordinate of the intersection points.

Parameters
[in]lineVecThe vector direction followed by the line
[in]offsetThe vector going from the sphere centre to a point on the line.
[in]sphereRadiusThe radius of the sphere.
Returns
pair of quadratic roots. Returns infinity if there is not intersection.

Definition at line 1070 of file meshprocessing.cpp.

◆ PrepareSnakeForCFD()

void PrepareSnakeForCFD ( const snake snakein,
double  distanceTol,
mesh meshgeom,
std::vector< double > &  holeCoords 
)

Prepares the snake to be used for CFD, removes duplicate points and triangulates it.

Parameters
[in]snakeinThe snakein
[in]distanceTolThe distance tolerance
meshgeomThe meshgeom
holeCoordsThe hole coordinates

Definition at line 191 of file meshprocessing.cpp.

◆ PseudoSurfaceAngle()

double PseudoSurfaceAngle ( const mesh meshin,
const std::array< int, 2 > &  surfInds 
)

Calculates the pseudo surface angle.

This pseudo angle is the dot product between the normal, i

Parameters
[in]meshinThe input mesh
[in]surfIndsThe surface indices
Returns
dot product between surface normals if facing outwards

Definition at line 567 of file meshprocessing.cpp.

◆ SurfaceCentroids()

std::vector<double> SurfaceCentroids ( const mesh meshin)

Generate a vector of coordinates of points at the surfaces pseudo centroid.

Parameters
[in]meshinThe input mesh
Returns
vector of coordinates with one coordinate inside each surface

Definition at line 908 of file meshprocessing.cpp.

◆ SurfaceInternalLayers()

std::vector<double> SurfaceInternalLayers ( const mesh meshin,
int  nLayers 
)

Returns points on edges between surface pseudo centroid and vertices.

Parameters
[in]meshinThe input mesh to process
[in]nLayersThe number of layers of points (0: only the centre, 1: layer surroundin the centre)
Exceptions
std::invalid_argumentif the number of layers is below 0.
Returns
Vector of points containing the centres and additional points.

Definition at line 940 of file meshprocessing.cpp.

◆ VolumeCentroids()

std::vector<double> VolumeCentroids ( const mesh meshin)

Generate a vector of coordinates of points at the volume pseudo centroid.

Parameters
[in]meshinThe input mesh
Returns
vector of coordinates with one coordinate inside each volume

Definition at line 837 of file meshprocessing.cpp.

◆ VolumeInternalLayers()

std::vector<double> VolumeInternalLayers ( const mesh meshin,
int  nLayers 
)

Returns points on edges between volume pseudo centroid and vertices.

Parameters
[in]meshinThe input mesh to process
[in]nLayersThe number of layers of points (0: only the centre, 1: layer surroundin the centre)
Exceptions
std::invalid_argumentif the number of layers is below 0.
Returns
Vector of points containing the centres and additional points.

Definition at line 869 of file meshprocessing.cpp.