Tri3#

class planestress.analysis.finite_element.Tri3(el_idx: int, el_tag: int, coords: npt.NDArray[np.float64], node_idxs: list[int], material: Material)[source]#

Bases: TriangularElement

Class for a three-noded linear triangular element.

Methods

b_matrix_jacobian

Calculates the B matrix and jacobian at an isoparametric point.

element_load_vector

Assembles the load vector for the element.

element_stiffness_matrix

Assembles the stiffness matrix for the element.

gauss_points

Gaussian weights and locations for n_point Gaussian integration.

get_element_results

Calculates various results for the finite element given nodal displacements.

get_triangulation

Returns a list of triangle indices for a Tri3 element.

iso_to_global

Converts a point in isoparametric coordinates to global coordinates.

nodal_isoparametric_coordinates

Returns the values of the isoparametric coordinates at the nodes.

shape_functions

Returns the shape functions at a point for a Tri3 element.

shape_functions_derivatives

Returns the derivatives of the shape functions at a point for a Tri3 element.

__init__(el_idx: int, el_tag: int, coords: npt.NDArray[np.float64], node_idxs: list[int], material: Material) None[source]#

Inits the Tri3 class.

Parameters:
  • el_idx (int) – Element index.

  • el_tag (int) – Element mesh tag.

  • coords (npt.NDArray[np.float64]) – A 2 x 3 numpy.ndarray of coordinates defining the element, i.e. [[x1, x2, x3], [y1, y2, y3]].

  • node_idxs (list[int]) – A list of node indexes defining the element, i.e. [idx1, idx2, idx3].

  • material (Material) – Material of the element.

static shape_functions(iso_coords: tuple[float, float, float]) ndarray[Any, dtype[float64]][source]#

Returns the shape functions at a point for a Tri3 element.

Parameters:

iso_coords (tuple[float, float, float]) – Location of the point in isoparametric coordinates.

Returns:

The values of the shape functions [N1, N2, N3].

Return type:

ndarray[Any, dtype[float64]]

static shape_functions_derivatives(iso_coords: tuple[float, float, float]) ndarray[Any, dtype[float64]][source]#

Returns the derivatives of the shape functions at a point for a Tri3 element.

Parameters:

iso_coords (tuple[float, float, float]) – Location of the point in isoparametric coordinates.

Returns:

The partial derivatives of the shape functions.

Return type:

ndarray[Any, dtype[float64]]

static nodal_isoparametric_coordinates() ndarray[Any, dtype[float64]][source]#

Returns the values of the isoparametric coordinates at the nodes.

Returns:

Values of the isoparametric coordinates at the nodes.

Return type:

ndarray[Any, dtype[float64]]

get_triangulation() list[tuple[int, int, int]][source]#

Returns a list of triangle indices for a Tri3 element.

Returns:

List of triangle indices.

Return type:

list[tuple[int, int, int]]

b_matrix_jacobian(iso_coords: tuple[float, float, float]) tuple[ndarray[Any, dtype[float64]], float]#

Calculates the B matrix and jacobian at an isoparametric point.

Parameters:

iso_coords (tuple[float, float, float]) – Location of the point in isoparametric coordinates.

Returns:

Derivatives of the shape function (B matrix) and value of the jacobian, (b_mat, j).

Return type:

tuple[ndarray[Any, dtype[float64]], float]

element_load_vector() ndarray[Any, dtype[float64]]#

Assembles the load vector for the element.

Returns:

Element load vector.

Return type:

ndarray[Any, dtype[float64]]

element_stiffness_matrix() ndarray[Any, dtype[float64]]#

Assembles the stiffness matrix for the element.

Returns:

Element stiffness matrix.

Return type:

ndarray[Any, dtype[float64]]

gauss_points(n_points: int) ndarray[Any, dtype[float64]]#

Gaussian weights and locations for n_point Gaussian integration.

Parameters:

n_points (int) – Number of gauss points.

Raises:

ValueError – If n_points is not 1, 3, 4 or 6.

Returns:

Gaussian weights and locations. For each gauss point - [weight, eta, xi, zeta].

Return type:

ndarray[Any, dtype[float64]]

get_element_results(u: ndarray[Any, dtype[float64]]) ElementResults#

Calculates various results for the finite element given nodal displacements.

Calculates the following:

  • Stresses at nodes

  • TODO

Parameters:

u (ndarray[Any, dtype[float64]]) – Displacement vector for the element.

Returns:

ElementResults object.

Return type:

ElementResults

iso_to_global(iso_coords: tuple[float, float, float]) tuple[float, float]#

Converts a point in isoparametric coordinates to global coordinates.

Parameters:

iso_coords (tuple[float, float, float]) – Location of the point in isoparametric coordinates.

Returns:

Location of the point in global coordinates (x, y).

Return type:

tuple[float, float]