WARPXM v1.10.0
Loading...
Searching...
No Matches
WmUnstructuredGeometry Class Reference

Class used for generating mesh geometry concerns. More...

#include <wmunstructuredgeometry.h>

Detailed Description

Class used for generating mesh geometry concerns.

Public Types

typedef std::pair< std::string, std::vector< int > > nodeGroupPair_t
 
typedef std::map< std::string, std::vector< int > > nodeGroupMap_t
 

Public Member Functions

 WmUnstructuredGeometry (const WmUnstructuredPatch &patch)
 Create unstructured geometry for a given patch.
 
 ~WmUnstructuredGeometry ()
 Destroy geometry.
 
void generateGeometry (real *coordinates, int *connectivity, int *neighborhood, int *orientations)
 Generates geometry concerns.
 
void addNodeGroup (const std::string name, const std::vector< int > &nodeIndexes)
 Add node group to geometry.
 
const realgetCoordinates () const
 Access coordinates.
 
const int * getConnectivity () const
 Access connectivity.
 
const int * getNeighborhood () const
 Access neighborhood.
 
const int * getOrientations () const
 Access orientations.
 
std::string getNodeGroupOfNodes (const std::vector< int > &nodes) const
 Identify which node group a set of nodes belong to.
 
int getFaceOrientation (const int elementIndex, const int faceIndex) const
 Get the orientation of a specific face with respect to neighbor.
 
real getVolume (const int elementIndex) const
 Get the volume of an element.
 
real getArea (const int elementIndex, const int faceIndex) const
 Get area of element face.
 
real getElementWidth (const int elementIndex) const
 Get rough estimate of element's characteristic width.
 
const int * getLocalNeighborhood (const int elementIndex) const
 Get list of element neighbors.
 
const std::type_info & getType () const
 Get typeid of floats used to store coordinates/geometry stuff.
 
int getPrimitiveID () const
 Get primitive ID for geometry.
 
int getFacePrimitiveID () const
 Get primitive's face primitive ID.
 
const realgetCentroid (const int elementIndex) const
 Get centroid of element.
 
const realgetNorm (const int elementIndex, const int faceIndex, const int normType) const
 Get norm of element.
 
int getNumElements () const
 Number of total elements in geometry.
 
int getNumNodes () const
 Number of total nodes in geometry.
 
int getNumDimsPerNode () const
 Number of dimensions used to store coordinates.
 
int getNumDimsPerPrimitive () const
 Number of dimensions used to express primtive.
 
int getNumNodesPerFace () const
 Number of nodes/vertexes on an element's face.
 
int getNumNodesPerElement () const
 Number of nodes/vertexes on an element.
 
int getNumFacesPerElement () const
 Number of faces on an element.
 
void getFaceNodes (const int elementIndex, const int faceIndex, int *faceNodes) const
 Get the nodes for a given element's face.
 
void getNodePositions (const int elementIndex, real *nodePositions) const
 Get the node/vertex positions for an element.
 
int getNumTransformCoordinates () const
 Get number of coordinate transform coordinates for this primitive.
 
void getTransformCoordinates (const int elementIndex, real tC[][3]) const
 Get the transform coordinates for an element.
 
void getPositionInElement (const int elementIndex, const real normalizedPosition[3], real position[3]) const
 Generate real/mesh space position from normalized position for given element.
 
int getJumpedElementIndex (const int elementIndex, const int neighborElementIndex, const int numJumps) const
 Get element index on other side of neighbor.
 
int getOppositeFace (const int faceIndex) const
 Get face index opposite of a given face index.
 
int getFaceIndex (const int elementIndex, const int neighborElementIndex) const
 Get the face index through which an neighbor element is found.
 
real generateVolumeJacobian (const int elementIndex, const real normalizedPosition[3]) const
 Generate the volume Jacobian (ie Jacobian or determinant of the Jacobian matrix given appropriate dimensions, I think) at a given point within an element.
 
void generateInverseJacobian (const int elementIndex, const real *normalizedPosition, real *inverseJacobian) const
 Generate the inverse Jacobian matrix for a given point within an element.
 
void getCentralPositionInElement (const int elementIndex, const real normalizedPosition[3], real position[3]) const
 Get the real/mesh position within an element w.r.t the element centroid.
 
void getFaceDistanceBetweenNeighborCentroids (const int fromElementIndex, const int toElementIndex, real dx[3]) const
 Get the distance between two elements through a face.
 
void get_face_centroid (const int element_index, const int face_index, real x[3]) const
 Get the centroid of element's face.
 
bool isEnclosedByElement (const int elementIndex, const real *x) const
 Check if the element encloses the given location.
 
int getEnclosingElementIndex (const real *x) const
 Find the element enclosing the given location.
 
void printElement (const int elementIndex) const
 Print element info to std::cout - used for DEBUG.
 
int num_boundary_faces () const
 Returns the number of boundary faces associated with the edge of a domain (i.e.
 
void reflect_centroid (const int element_index, const int face_index, real reflected_centroid[3]) const
 

Static Public Member Functions

static int getNumTransformCoordinates (const int primitiveID)
 Number of coordinate transform coordinates.
 
static void getTransformCoordinates (const int primitiveID, const real *nodePositions, real tC[])
 Convert node/vertex positions into coordinate transform coordinates.
 

Protected Attributes

const WmUnstructuredPatch_patch
 Patch associated with geometry.
 
std::string _primitive
 Name of element primitive.
 
int _primitiveID
 PrimitiveID.
 
int _numElements
 Number of total elements in patch.
 
int _numNodes
 Number of total nodes in patch.
 
int * _connectivity
 Connectivity of elements [_numElements X _numNodesPerElement].
 
int * _neighborhood
 Neighborhood of elements [_numElements X _numFacesPerElement].
 
int * _orientations
 Orientations of faces [_numElements X _numFacesPerElement].
 
real_coordinates
 Node/vertex positions [_numNodes X _numDimsPerNode].
 
real_transformCoordinates
 Coordinate transform coordinates [_numElements X _numTransformCoordinates X 3].
 
real_centroids
 Centroids of elements [_numElements X 3].
 
real_norms
 Face normals for elements [_numElements X 3 X 3].
 
real_volumes
 Volumes of elements [_numElements].
 
real_areas
 Areas of element faces [_numElements X _numFacesPerElement].
 
nodeGroupMap_t _nodeGroups
 Node groups associated with geometry/patch.
 
int _num_boundary_faces
 Number of boundary faces between patch and unknown where neighbor = -1.
 
int _numDimsPerNode
 Number of dimensions used to store nodes/vertexes.
 
int _numDimsPerPrimitive
 Number of dimensions for a given primitive.
 
int _numTransformCoordinates
 Number of coordinate transform coordinates, usually the same as _numNodesPerElement.
 
int _numNodesPerElement
 Number of nodes/vertexes for a given primitive.
 
int _numFacesPerElement
 Number of faces on an element.
 
int _numNodesPerFace
 Number of element nodes/vertexes that make up each face.
 

Member Typedef Documentation

◆ nodeGroupMap_t

typedef std::map<std::string, std::vector<int> > WmUnstructuredGeometry::nodeGroupMap_t

◆ nodeGroupPair_t

typedef std::pair<std::string, std::vector<int> > WmUnstructuredGeometry::nodeGroupPair_t

Constructor & Destructor Documentation

◆ WmUnstructuredGeometry()

WmUnstructuredGeometry::WmUnstructuredGeometry ( const WmUnstructuredPatch patch)

Create unstructured geometry for a given patch.

◆ ~WmUnstructuredGeometry()

WmUnstructuredGeometry::~WmUnstructuredGeometry ( )

Destroy geometry.

Member Function Documentation

◆ addNodeGroup()

void WmUnstructuredGeometry::addNodeGroup ( const std::string  name,
const std::vector< int > &  nodeIndexes 
)

Add node group to geometry.

Node groups are used to identify faces (e.g. for use with boundary conditions)

Note: node group will be sorted in geometry class

Parameters
nameNode group name
nodeIndexesList of nodes that make up node group

◆ generateGeometry()

void WmUnstructuredGeometry::generateGeometry ( real coordinates,
int *  connectivity,
int *  neighborhood,
int *  orientations 
)

Generates geometry concerns.

Includes: coordinate transforms, face normals, volumes, face areas and element centroids

Parameters
coordinatesCoordinates for mesh nodes/vertexes
connecitivityConnectivity array for mesh elements
neighborhoodNeighborhood array interlinking elements
orientationsOrienatation array describing inter-element orientations

◆ generateInverseJacobian()

void WmUnstructuredGeometry::generateInverseJacobian ( const int  elementIndex,
const real normalizedPosition,
real inverseJacobian 
) const

Generate the inverse Jacobian matrix for a given point within an element.

Parameters
elementIndexElement of interest
normalizedPositionArray [3] containing normalized position within element
inverseJacobianArray [9] to store inverse Jacobian matrix

◆ generateVolumeJacobian()

real WmUnstructuredGeometry::generateVolumeJacobian ( const int  elementIndex,
const real  normalizedPosition[3] 
) const

Generate the volume Jacobian (ie Jacobian or determinant of the Jacobian matrix given appropriate dimensions, I think) at a given point within an element.

Parameters
elementIndexElement of interest
normalizedPositionArray [3] containing normalized position within element
Returns
Volume Jacobian at point within element

◆ get_face_centroid()

void WmUnstructuredGeometry::get_face_centroid ( const int  element_index,
const int  face_index,
real  x[3] 
) const

Get the centroid of element's face.

Parameters
element_index
face_index
xArray [3] to store position of face centroid

◆ getArea()

real WmUnstructuredGeometry::getArea ( const int  elementIndex,
const int  faceIndex 
) const
inline

Get area of element face.

Parameters
elementIndexElement of interest
faceIndexElement face of interest
Returns
Area of element face

◆ getCentralPositionInElement()

void WmUnstructuredGeometry::getCentralPositionInElement ( const int  elementIndex,
const real  normalizedPosition[3],
real  position[3] 
) const

Get the real/mesh position within an element w.r.t the element centroid.

Parameters
elementIndexElement of interest
normalizedPositionArray [3] containing normalized position within element
positionArray [3] to store central position

◆ getCentroid()

const real * WmUnstructuredGeometry::getCentroid ( const int  elementIndex) const
inline

Get centroid of element.

Warning: no range checking!

Parameters
elementIndexElement of interest
Returns
Array [3] containing element centroid

◆ getConnectivity()

const int * WmUnstructuredGeometry::getConnectivity ( ) const
inline

Access connectivity.

Returns
Connectivity of elements [numElements X numNodesPerElement]

◆ getCoordinates()

const real * WmUnstructuredGeometry::getCoordinates ( ) const
inline

Access coordinates.

Returns
Coordinates of mesh nodes/vertexes [numNodes X numDimsPerNode]

◆ getElementWidth()

real WmUnstructuredGeometry::getElementWidth ( const int  elementIndex) const

Get rough estimate of element's characteristic width.

Parameters
elementIndexElement of interest
Returns
'Width' of element

◆ getEnclosingElementIndex()

int WmUnstructuredGeometry::getEnclosingElementIndex ( const real x) const

Find the element enclosing the given location.

Note that this call is going to run through all of the elements and compare each one Therefore this is a very slow function and should be called only a few times at most per sim

Note that this function only works if 1D) Line is along x 2D) Triangle is in xy plane 3D) Not yet setup

Parameters
xArray[3] of location
Returns
Element Index of element enclosing x (-1 if no element was found)

◆ getFaceDistanceBetweenNeighborCentroids()

void WmUnstructuredGeometry::getFaceDistanceBetweenNeighborCentroids ( const int  fromElementIndex,
const int  toElementIndex,
real  dx[3] 
) const

Get the distance between two elements through a face.

Note: This call is required for periodic boundary conditions where two elements are on opposite sides of a mesh but require 'dx' for reconstruction purposes Note: This is a slow call compared to subracting one centroid from another, be sure you need it before using it

dx = x_toIndex - x_fromIndex

Parameters
fromElementIndexElement to start from
toElementIndexElement to end at
dxArray [3] to store dx between elements through face

◆ getFaceIndex()

int WmUnstructuredGeometry::getFaceIndex ( const int  elementIndex,
const int  neighborElementIndex 
) const

Get the face index through which an neighbor element is found.

Parameters
elementIndexElement of interest
neighborElementIndexNeighbor element on other side of face
Returns
Face index separating elements

◆ getFaceNodes()

void WmUnstructuredGeometry::getFaceNodes ( const int  elementIndex,
const int  faceIndex,
int *  faceNodes 
) const

Get the nodes for a given element's face.

Parameters
elementIndexElement of interest
faceIndexElement face of interest
faceNodesArray [_numNodesPerFace] to store node/vertex indexes

◆ getFaceOrientation()

int WmUnstructuredGeometry::getFaceOrientation ( const int  elementIndex,
const int  faceIndex 
) const
inline

Get the orientation of a specific face with respect to neighbor.

Parameters
elementIndexElement of interest
faceIndexFace index through which neighbor is oriented
Returns
Orientation of neighbor element

◆ getFacePrimitiveID()

int WmUnstructuredGeometry::getFacePrimitiveID ( ) const

Get primitive's face primitive ID.

Will return one of WMPRIMITIVE_LINE, WMPRIMITIVE_TRIANGLE, WMPRIMITIVE_QUADRILATERAL, WMPRIMITIVE_TETRAHEDRON, WMPRIMITIVE_HEXAHEDRON, or WMPRIMITIVE_POINT

Returns
Primitive ID

◆ getJumpedElementIndex()

int WmUnstructuredGeometry::getJumpedElementIndex ( const int  elementIndex,
const int  neighborElementIndex,
const int  numJumps 
) const

Get element index on other side of neighbor.

I think this is used for moving through ghost cells in boundary conditions, but is... idiotic...

Note: Recursive call Note: Only works for lines/quad/hexahedrons

TODO: This function needs to be removed/phased out

Parameters
elementIndexElement of interest
neighborElementIndexNeighbor element to jump over
numJumpsNumber of times to jump across neighbor
Returns
Element index numJumps beyond neighborElementIndex

◆ getLocalNeighborhood()

const int * WmUnstructuredGeometry::getLocalNeighborhood ( const int  elementIndex) const
inline

Get list of element neighbors.

Warning: no range checking!

Parameters
elementIndexElement of interest
Returns
Array of neighbor elements

◆ getNeighborhood()

const int * WmUnstructuredGeometry::getNeighborhood ( ) const
inline

Access neighborhood.

Returns
Neighborhood of elements [numElements X numFacesPerElement]

◆ getNodeGroupOfNodes()

std::string WmUnstructuredGeometry::getNodeGroupOfNodes ( const std::vector< int > &  nodes) const

Identify which node group a set of nodes belong to.

Used for identifying faces contained by node group

Parameters
nodesList of nodes
Returns
Name of node group containing nodes, returns "" if nodes do not belong to node group

◆ getNodePositions()

void WmUnstructuredGeometry::getNodePositions ( const int  elementIndex,
real nodePositions 
) const

Get the node/vertex positions for an element.

Parameters
elementIndexElement of interest
nodePositionsArray [3 X numNodesPerElement] to store node/vertex positions

◆ getNorm()

const real * WmUnstructuredGeometry::getNorm ( const int  elementIndex,
const int  faceIndex,
const int  normType 
) const
inline

Get norm of element.

Warning: no range checking!

norm type 0: face normal norm type 1: first transverse vector norm type 2: second transverse vector

Parameters
elementIndexElement of interest
faceIndexElement of interest
normTypeNorm type of interese (see above)
Returns
Array [3] containing face norm of given type

◆ getNumDimsPerNode()

int WmUnstructuredGeometry::getNumDimsPerNode ( ) const
inline

Number of dimensions used to store coordinates.

Returns
Number of dims used to store coordinates

◆ getNumDimsPerPrimitive()

int WmUnstructuredGeometry::getNumDimsPerPrimitive ( ) const
inline

Number of dimensions used to express primtive.

1D: Line = 1 2D: Triangle, Quad = 2 3D: Tetrahedron, Hexahedron = 3

Returns
Number of dims of primitive

◆ getNumElements()

int WmUnstructuredGeometry::getNumElements ( ) const
inline

Number of total elements in geometry.

Returns
Number of total elements

◆ getNumFacesPerElement()

int WmUnstructuredGeometry::getNumFacesPerElement ( ) const
inline

Number of faces on an element.

Returns
Number of faces on an element

◆ getNumNodes()

int WmUnstructuredGeometry::getNumNodes ( ) const
inline

Number of total nodes in geometry.

Returns
Number of total nodes

◆ getNumNodesPerElement()

int WmUnstructuredGeometry::getNumNodesPerElement ( ) const
inline

Number of nodes/vertexes on an element.

Returns
Number of nodes/vertexes on an element

◆ getNumNodesPerFace()

int WmUnstructuredGeometry::getNumNodesPerFace ( ) const
inline

Number of nodes/vertexes on an element's face.

Returns
Number of nodes/vertexes on an element's face

◆ getNumTransformCoordinates() [1/2]

int WmUnstructuredGeometry::getNumTransformCoordinates ( ) const
inline

Get number of coordinate transform coordinates for this primitive.

Note: Coordinate transform converts isoparametric primitive (normalized space) to geometric primitive (real/mesh space)

Returns
Number of coordinate transform coordinates

◆ getNumTransformCoordinates() [2/2]

static int WmUnstructuredGeometry::getNumTransformCoordinates ( const int  primitiveID)
static

Number of coordinate transform coordinates.

Note: Coordinate transform converts isoparametric primitive (normalized space) to geometric primitive (real/mesh space)

Note: Static call

Parameters
primitiveIDPrimitive of interest
Returns
Number of coordinate transform coordinates

◆ getOppositeFace()

int WmUnstructuredGeometry::getOppositeFace ( const int  faceIndex) const

Get face index opposite of a given face index.

Another ghost cell function that shouldn't exist

Note: Static call Note: Only works for lines/quad/hexahedrons

TODO: This function needs to be removed/phased out

Parameters
faceIndexFace index
neighborElementIndexNeighbor element to jump over
numJumpsNumber of times to jump across neighbor
Returns
Opposite face index

◆ getOrientations()

const int * WmUnstructuredGeometry::getOrientations ( ) const
inline

Access orientations.

Returns
Orientation of elements [numElements X numFacesPerElement]

◆ getPositionInElement()

void WmUnstructuredGeometry::getPositionInElement ( const int  elementIndex,
const real  normalizedPosition[3],
real  position[3] 
) const

Generate real/mesh space position from normalized position for given element.

Parameters
elementIndexElement of interest
normalizedPositionArray [3] of position in isoparametric space (normalized space)
positionArray [3] of position in real/mesh space

◆ getPrimitiveID()

int WmUnstructuredGeometry::getPrimitiveID ( ) const
inline

Get primitive ID for geometry.

Will return one of WMPRIMITIVE_LINE, WMPRIMITIVE_TRIANGLE, WMPRIMITIVE_QUADRILATERAL, WMPRIMITIVE_TETRAHEDRON, WMPRIMITIVE_HEXAHEDRON, or WMPRIMITIVE_POINT

Returns
Primitive ID

◆ getTransformCoordinates() [1/2]

void WmUnstructuredGeometry::getTransformCoordinates ( const int  elementIndex,
real  tC[][3] 
) const

Get the transform coordinates for an element.

Parameters
elementIndexElement of interest
tCArray [numTransformCoordinates][3] to store transformCoordinates

◆ getTransformCoordinates() [2/2]

static void WmUnstructuredGeometry::getTransformCoordinates ( const int  primitiveID,
const real nodePositions,
real  tC[] 
)
static

Convert node/vertex positions into coordinate transform coordinates.

Note: Coordinate transform converts isoparametric primitive (normalized space) to geometric primitive (real/mesh space)

Note: Static call

Parameters
primitiveIDPrimitive of interest
nodePositionsNode positions making up element
tCArray [_numTransformCoordinates][3] containing coordinate transform coordinates for element

◆ getType()

const std::type_info & WmUnstructuredGeometry::getType ( ) const
inline

Get typeid of floats used to store coordinates/geometry stuff.

Returns
Typeid of floats

◆ getVolume()

real WmUnstructuredGeometry::getVolume ( const int  elementIndex) const
inline

Get the volume of an element.

Parameters
elementIndexElement of interest
Returns
Volume of element

◆ isEnclosedByElement()

bool WmUnstructuredGeometry::isEnclosedByElement ( const int  elementIndex,
const real x 
) const

Check if the element encloses the given location.

Parameters
elementIndexElement to test
xArray[3] of location
Returns
Does the element enclose the given point?

◆ num_boundary_faces()

int WmUnstructuredGeometry::num_boundary_faces ( ) const
inline

Returns the number of boundary faces associated with the edge of a domain (i.e.

neighbor = -1)

Returns
Number of boundary faces for patch with all attached node groups

◆ printElement()

void WmUnstructuredGeometry::printElement ( const int  elementIndex) const

Print element info to std::cout - used for DEBUG.

Parameters
elementIndexElement of interest

◆ reflect_centroid()

void WmUnstructuredGeometry::reflect_centroid ( const int  element_index,
const int  face_index,
real  reflected_centroid[3] 
) const

Member Data Documentation

◆ _areas

real* WmUnstructuredGeometry::_areas
protected

Areas of element faces [_numElements X _numFacesPerElement].

◆ _centroids

real* WmUnstructuredGeometry::_centroids
protected

Centroids of elements [_numElements X 3].

◆ _connectivity

int* WmUnstructuredGeometry::_connectivity
protected

Connectivity of elements [_numElements X _numNodesPerElement].

◆ _coordinates

real* WmUnstructuredGeometry::_coordinates
protected

Node/vertex positions [_numNodes X _numDimsPerNode].

◆ _neighborhood

int* WmUnstructuredGeometry::_neighborhood
protected

Neighborhood of elements [_numElements X _numFacesPerElement].

◆ _nodeGroups

nodeGroupMap_t WmUnstructuredGeometry::_nodeGroups
protected

Node groups associated with geometry/patch.

◆ _norms

real* WmUnstructuredGeometry::_norms
protected

Face normals for elements [_numElements X 3 X 3].

◆ _num_boundary_faces

int WmUnstructuredGeometry::_num_boundary_faces
protected

Number of boundary faces between patch and unknown where neighbor = -1.

◆ _numDimsPerNode

int WmUnstructuredGeometry::_numDimsPerNode
protected

Number of dimensions used to store nodes/vertexes.

◆ _numDimsPerPrimitive

int WmUnstructuredGeometry::_numDimsPerPrimitive
protected

Number of dimensions for a given primitive.

◆ _numElements

int WmUnstructuredGeometry::_numElements
protected

Number of total elements in patch.

◆ _numFacesPerElement

int WmUnstructuredGeometry::_numFacesPerElement
protected

Number of faces on an element.

◆ _numNodes

int WmUnstructuredGeometry::_numNodes
protected

Number of total nodes in patch.

◆ _numNodesPerElement

int WmUnstructuredGeometry::_numNodesPerElement
protected

Number of nodes/vertexes for a given primitive.

◆ _numNodesPerFace

int WmUnstructuredGeometry::_numNodesPerFace
protected

Number of element nodes/vertexes that make up each face.

◆ _numTransformCoordinates

int WmUnstructuredGeometry::_numTransformCoordinates
protected

Number of coordinate transform coordinates, usually the same as _numNodesPerElement.

◆ _orientations

int* WmUnstructuredGeometry::_orientations
protected

Orientations of faces [_numElements X _numFacesPerElement].

◆ _patch

const WmUnstructuredPatch& WmUnstructuredGeometry::_patch
protected

Patch associated with geometry.

◆ _primitive

std::string WmUnstructuredGeometry::_primitive
protected

Name of element primitive.

◆ _primitiveID

int WmUnstructuredGeometry::_primitiveID
protected

PrimitiveID.

Is either: WMPRIMITIVE_LINE, WMPRIMITIVE_TRIANGLE, WMPRIMITIVE_QUADRILATERAL, WMPRIMITIVE_TETRAHEDRON, WMPRIMITIVE_HEXAHEDRON, or WMPRIMITIVE_POINT

◆ _transformCoordinates

real* WmUnstructuredGeometry::_transformCoordinates
protected

Coordinate transform coordinates [_numElements X _numTransformCoordinates X 3].

◆ _volumes

real* WmUnstructuredGeometry::_volumes
protected

Volumes of elements [_numElements].


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