spharapy.trimesh: Triangular Mesh Data

Triangular mesh data

This module provides a class for storing triangular meshes. Attributes of the triangular mesh can be determined. In addition, methodes are available to derive further information from the triangular grid.

class spharapy.trimesh.TriMesh(trilist, vertlist)[source]

Bases: object

Triangular mesh class

This class can be used to store data to define a triangular mesh and it provides atributes and methodes to derive further information about the triangular mesh.

Parameters
trilist: array, shape (n_triangles, 3)

List of triangles, each row of the array contains the edges of a triangle. The edges of the triangles are defined by the indices to the list of vertices. The index of the first vertex is 0. The number of triangles is n_triangles.

vertlist: array, shape (n_points, 3)

List of coordinates x, y, z which describes the positions of the vertices.

Attributes
trilist: array, shape (n_triangles, 3)

List of triangles of the mesh.

vertlist: array, shape (n_points, 3)

List of coordinates of the vertices

adjacent_tri(vertex_index=0)[source]

All triangles with the given vertex

The method determined all triangles of the triangular mesh that contain the given vertex.

Parameters
vertex_indexinteger

Index of the vertex for which the adjacent vertices are to be determined. The index must be in the range 0 to number of vertices - 1.

Returns
tri_with_vertexarray, shape (3, n)

List of triangles containing the given vertex.

is_edge(vertex1_index, vertex2_index)[source]

Are 2 vertices connected by an edge

The method determines whether two vertices are connected by an edge in the triangle mesh and if so, whether it is an internal edge or a boundary edge.

Parameters
vertex1_index, vertex2_indexinteger

Indeces of the two vertices. The index must be in the range 0 to number of vertices - 1.

Returns
is_edgeinteger

0 if vertex1 and vertex2 are not connected by a single edge, 1 if vertex1 and vertex2 are connected by a boundary edge, 2 if vertex1 and vertex2 are connected by an internal edge.

laplacianmatrix(mode='inv_euclidean')[source]

Compute a laplacian matrix for a triangular mesh

The method creates a laplacian matrix for a triangular mesh using different weighting function.

Parameters
mode{‘unit’, ‘inv_euclidean’, ‘half_cotangent’}, optional

The methods for determining the edge weights. Using the option ‘unit’ all edges of the mesh are weighted by unit weighting function, the result is an adjacency matrix. The option ‘inv_euclidean’ results in edge weights corresponding to the inverse Euclidean distance of the edge lengths. The option ‘half_cotangent’ uses the half of the cotangent of the two angles opposed to an edge as weighting function. the default weighting function is ‘inv_euclidean’.

Returns
laplacianmatrixarray, shape (n_points, n_points)

Matrix, which contains the discrete laplace operator for data defined at the vertices of a triangular mesh. The number of vertices of the triangular mesh is n_points.

Examples

>>> from spharapy import trimesh as tm
>>> testtrimesh = tm.TriMesh([[0, 1, 2]], [[1., 0., 0.], [0., 2., 0.],
...                                        [0., 0., 3.]])
>>> testtrimesh.laplacianmatrix(mode='inv_euclidean')
array([[ 0.76344136, -0.4472136 , -0.31622777],
       [-0.4472136 ,  0.72456369, -0.2773501 ],
       [-0.31622777, -0.2773501 ,  0.59357786]])
massmatrix(mode='normal')[source]

Mass matrix of a triangular mesh

The method determines a mass matrix of a triangular mesh.

Parameters
mode{‘normal’, ‘lumped’}, optional

The mode parameter can be used to select whether a normal mass matrix or a lumped mass matrix is to be determined.

Returns
massmatrixarray, shape (n_points, n_points)

Symmetric matrix, which contains the mass values for each edge and vertex for the FEM approch. The number of vertices of the triangular mesh is n_points.

References

[VL07][DZMC07][ZvKD07]

Examples

>>> from spharapy import trimesh as tm
>>> testtrimesh = tm.TriMesh([[0, 1, 2]], [[1., 0., 0.], [0., 2., 0.],
...                                        [0., 0., 3.]])
>>> testtrimesh.massmatrix()
array([[ 0.58333333,  0.29166667,  0.29166667],
       [ 0.29166667,  0.58333333,  0.29166667],
       [ 0.29166667,  0.29166667,  0.58333333]])
one_ring_neighborhood(vertex_index=0)[source]

The 1 ring neighborhood of a vertex

The method determines all adjacent vertices of a vertex that is given by its index, the so called 1 ring neighborhood.

Parameters
vertex_indexinteger

Index of the vertex for which the adjacent vertices are to be determined. The index must be in the range 0 to number of vertices - 1.

Returns
one_ring_neighborhoodarray, shape (1, n)

Array of indexes on vertices adjacent to a given vertex.

remove_vertices(vertex_index_list)[source]

Remove vertices from a triangular mesh

The method removes vertices from a triangle mesh. The Half-edge Collapse method is used. The positions of the remaining vertices are not affected and they are retriangulated.

Parameters
vertex_index_listvector of ints

Indices of the vertices to remove from the mesh. The indices must be in the range 0 to number of vertices - 1.

Returns
triangsamplestrimesh object

A trimesh object from the package spharapy, where the given vertices are removed.

stiffnessmatrix()[source]

Stiffness matrix of a triangular mesh

The method determines a stiffness matrix of a triangular mesh.

Returns
stiffmatrixarray, shape (n_points, n_points)

Symmetric matrix, which contains the stiffness values for each edge and vertex for the FEM approch. The number of vertices of the triangular mesh is n_points.

References

[VL07]

Examples

>>> from spharapy import trimesh as tm
>>> testtrimesh = tm.TriMesh([[0, 1, 2]], [[1., 0., 0.], [0., 2., 0.],
...                                        [0., 0., 3.]])
>>> testtrimesh.stiffnessmatrix()
array([[-0.92857143,  0.64285714,  0.28571429],
       [ 0.64285714, -0.71428571,  0.07142857],
       [ 0.28571429,  0.07142857, -0.35714286]])
property trilist

Get or set the list of triangles.

Setting the list of triangles will simultaneously check if the triangle list is in the correct format.

property vertlist

Get or set the list of vertices.

Setting the list of triangles will simultaneously check if the vertice list is in the correct format.

weightmatrix(mode='inv_euclidean')[source]

Compute a weight matrix for a triangular mesh

The method creates a weighting matrix for the edges of a triangular mesh using different weighting function.

Parameters
mode{‘unit’, ‘inv_euclidean’, ‘half_cotangent’}, optional

The parameter mode specifies the method for determining the edge weights. Using the option ‘unit’ all edges of the mesh are weighted by unit weighting function, the result is an adjacency matrix. The option ‘inv_euclidean’ results in edge weights corresponding to the inverse Euclidean distance of the edge lengths. The option ‘half_cotangent’ uses the half of the cotangent of the two angles opposed to an edge as weighting function. the default weighting function is ‘inv_euclidean’.

Returns
weightmatrixarray, shape (n_points, n_points)

Symmetric matrix, which contains the weight of the edges between adjacent vertices. The number of vertices of the triangular mesh is n_points.

Examples

>>> from spharapy import trimesh as tm
>>> testtrimesh = tm.TriMesh([[0, 1, 2]], [[1., 0., 0.], [0., 2., 0.],
...                                        [0., 0., 3.]])
>>> testtrimesh.weightmatrix(mode='inv_euclidean')
array([[ 0.        ,  0.4472136 ,  0.31622777],
       [ 0.4472136 ,  0.        ,  0.2773501 ],
       [ 0.31622777,  0.2773501 ,  0.        ]])
spharapy.trimesh.angles_triangle(vertex1, vertex2, vertex3)[source]

Estimate the three internal angles of a triangle given by three vertices

Parameters
vertex1array, shape (1, 3)
vertex2array, shape (1, 3)
vertex3array, shape (1, 3)
Returns
anglesarray, shape (1, 3)

Internal angles of the triangle given by the three vertices.

Examples

>>> from spharapy import trimesh as tm
>>> tm.angles_triangle([1, 0, 0], [0, 1, 0], [0, 0, 1])
array([1.04719755, 1.04719755, 1.04719755])
spharapy.trimesh.area_triangle(vertex1, vertex2, vertex3)[source]

Estimate the area of a triangle given by three vertices

The area of the triangle given by three vertices is calculated by the half cross product formula.

Parameters
vertex1array, shape (1, 3)
vertex2array, shape (1, 3)
vertex3array, shape (1, 3)
Returns
triangleareafloat

Area of the triangle given by the three vertices.

Examples

>>> from spharapy import trimesh as tm
>>> tm.area_triangle([1, 0, 0], [0, 1, 0], [0, 0, 1])
0.8660254037844386
spharapy.trimesh.side_lens_triangle(vertex1, vertex2, vertex3)[source]

Estimate the three side length of a triangle given by three vertices

Parameters
vertex1array, shape (1, 3)
vertex2array, shape (1, 3)
vertex3array, shape (1, 3)
Returns
side_lensarray, shape (1, 3)

Side lengths of the triangle given by the three vertices.

Examples

>>> from spharapy import trimesh as tm
>>> tm.side_lens_triangle([1, 0, 0], [0, 1, 0], [0, 0, 1])
array([1.41421356, 1.41421356, 1.41421356])