LoopStructural.interpolators.P2UnstructuredTetMesh#

class LoopStructural.interpolators.P2UnstructuredTetMesh(nodes: ndarray, elements: ndarray, neighbours: ndarray, aabb_nsteps=None)#

Bases: UnStructuredTetMesh

An unstructured mesh defined by nodes, elements and neighbours An axis aligned bounding box (AABB) is used to speed up finding which tetra a point is in. The aabb grid is calculated so that there are approximately 10 tetra per element.

Parameters:
  • nodes (array or array like) – container of vertex locations

  • elements (array or array like, dtype cast to long) – container of tetra indicies

  • neighbours (array or array like, dtype cast to long) – array containing element neighbours

  • aabb_nsteps (list, optional) – force nsteps for aabb, by default None

__init__(nodes: ndarray, elements: ndarray, neighbours: ndarray, aabb_nsteps=None)#

An unstructured mesh defined by nodes, elements and neighbours An axis aligned bounding box (AABB) is used to speed up finding which tetra a point is in. The aabb grid is calculated so that there are approximately 10 tetra per element.

Parameters:
  • nodes (array or array like) – container of vertex locations

  • elements (array or array like, dtype cast to long) – container of tetra indicies

  • neighbours (array or array like, dtype cast to long) – array containing element neighbours

  • aabb_nsteps (list, optional) – force nsteps for aabb, by default None

Methods

__init__(nodes, elements, neighbours[, ...])

An unstructured mesh defined by nodes, elements and neighbours An axis aligned bounding box (AABB) is used to speed up finding which tetra a point is in.

evaluate_d2(pos, prop)

Evaluate the second derivative of the interpolant d2x, dxdy, d2y, dxdz dydz d2dz

evaluate_gradient(pos, property_array)

Evaluate the gradient of an interpolant at the locations

evaluate_shape(locations)

Convenience function returning barycentric coords

evaluate_shape_d2(indexes)

evaluate second derivatives of shape functions in s and t

evaluate_shape_derivatives(locations[, elements])

compute dN/ds (1st row), dN/dt(2nd row)

evaluate_value(pos, property_array)

Evaluate value of interpolant

get_element_for_location(points)

Determine the tetrahedron from a numpy array of points

get_element_gradient_for_location(pos)

Get the gradient of the tetra for a location

get_element_gradients([elements])

Get the gradients of all tetras

get_elements()

get_neighbours()

This function goes through all of the elements in the mesh and assembles a numpy array with the neighbours for each element

get_quadrature_points([npts])

Calculate the quadrature points for the triangle using 3 points these points are at the barycentric coordinates of (1/6,1/6), (1/6,2/3), (2/3,1/6) All points are weighted equally at 1/6

inside(pos)

Check if a position is inside the support

onGeometryChange()

Called when the geometry changes

Attributes

barycentre

Return the number of dimensions

dimension

element_size

Calculate the volume of a tetrahedron using the 4 corners volume = abs(det(A))/6 where A is the jacobian of the corners

elements

Return the elements

n_cells

n_elements

Return the number of elements

n_nodes

Return the number of points

nodes

Return the nodes

ntetra

shared_element_norm

Get the normal to all of the shared elements

shared_element_size

Get the area of the share triangle

vtk

property barycentre#

Return the number of dimensions

property element_size#

Calculate the volume of a tetrahedron using the 4 corners volume = abs(det(A))/6 where A is the jacobian of the corners

Returns:

_type_ – _description_

property elements#

Return the elements

evaluate_d2(pos: ndarray, prop: ndarray) ndarray#

Evaluate the second derivative of the interpolant d2x, dxdy, d2y, dxdz dydz d2dz

Parameters:
  • array (prop - numpy) – locations

  • array – property values at nodes

evaluate_gradient(pos: ndarray, property_array: ndarray) ndarray#

Evaluate the gradient of an interpolant at the locations

Parameters:
  • array (pos - numpy) – locations

  • string (prop -) – property to evaluate

evaluate_shape(locations: ndarray)#

Convenience function returning barycentric coords

evaluate_shape_d2(indexes: ndarray) ndarray#

evaluate second derivatives of shape functions in s and t

Parameters:

indexes (np.ndarray) – array of indexes

Returns:

np.array – array of second derivative shape function

evaluate_shape_derivatives(locations: ndarray, elements: ndarray | None = None) ndarray#

compute dN/ds (1st row), dN/dt(2nd row)

Parameters:
  • locations (np.array) – location (n,3) array

  • elements (np.array, optional) – indexes to calculate shape function for. Used when evaluating quad points on faces as two tetra hold the point by default None When it is none, the index is calculated from the location

Returns:

np.array – array of shape paramters

evaluate_value(pos: ndarray, property_array: ndarray) ndarray#

Evaluate value of interpolant

Parameters:
  • array (pos - numpy) – locations

  • string (prop -) – property name

get_element_for_location(points: ndarray) Tuple#

Determine the tetrahedron from a numpy array of points

Parameters:

pos (np.array) –

get_element_gradient_for_location(pos)#

Get the gradient of the tetra for a location

Parameters:

pos

get_element_gradients(elements=None)#

Get the gradients of all tetras

Parameters:

elements

get_neighbours()#

This function goes through all of the elements in the mesh and assembles a numpy array with the neighbours for each element

get_quadrature_points(npts: int = 3)#

Calculate the quadrature points for the triangle using 3 points these points are at the barycentric coordinates of (1/6,1/6), (1/6,2/3), (2/3,1/6) All points are weighted equally at 1/6

Parameters:

npts (int, optional) – _description_, by default 3

Returns:

_type_ – _description_

inside(pos)#

Check if a position is inside the support

property n_elements#

Return the number of elements

property n_nodes#

Return the number of points

property nodes#

Return the nodes

onGeometryChange()#

Called when the geometry changes

property shared_element_norm#

Get the normal to all of the shared elements

property shared_element_size#

Get the area of the share triangle