pysdic.compute_shape_functions#
- compute_shape_functions(natural_coordinates, element_type, return_derivatives=False, *, default=0.0)[source]#
Compute the shape functions for a given
element_typeat specifiednatural_coordinates.Note
Input
natural_coordinateswill be converted tonumpy.float64.Output arrays will be
numpy.float64.
- Parameters:
natural_coordinates (ArrayLike) – Natural coordinates where to evaluate the shape functions. The array must have shape \((N_{p}, K)\), where \(N_{p}\) is the number of points to evaluate and \(K\) is the topological dimension of the element.
element_type (
str) –Type of the finite element. Supported types are:
'segment_2': 2-node segment element.'segment_3': 3-node segment element.'triangle_3': 3-node triangle element.'triangle_6': 6-node triangle element.'quadrangle_4': 4-node quadrangle element.'quadrangle_8': 8-node quadrangle element.
return_derivatives (
bool, optional) – IfTrue, the function will also return the first derivatives of the shape functions with respect to the natural coordinates. By default,False.default (Real, optional) – Default value to assign to shape functions when the input natural coordinates are out of the valid range (i.e., not in the valid range for the specified element type). By default,
0.0.
- Returns:
shape_functions (
numpy.ndarray) – Shape functions evaluated at the given natural coordinates. The returned array has shape \((N_{p}, N_{vpe})\), where each row corresponds to a point and each column to a node, and \(N_{vpe}\) is the number of vertices per element for the specified element type.shape_function_derivatives (
numpy.ndarray, optional) – Ifreturn_derivativesisTrue, the function also returns an array of the first derivatives of the shape functions with respect to the natural coordinates. The returned array has shape \((N_{p}, N_{vpe}, K)\).
- Return type:
See also
pysdic.segment_2_shape_functionsShape functions for 2-node segment 1D-elements.
pysdic.compute_segment_3_shape_functionsShape functions for 3-node segment 1D-elements.
pysdic.compute_triangle_3_shape_functionsShape functions for 3-node triangle 2D-elements.
pysdic.compute_triangle_6_shape_functionsShape functions for 6-node triangle 2D-elements.
pysdic.compute_quadrangle_4_shape_functionsShape functions for 4-node quadrangle 2D-elements.
pysdic.compute_quadrangle_8_shape_functionsShape functions for 8-node quadrangle 2D-elements.
Examples
Compute shape functions for a 6-node triangle element for 3 valid points and 1 invalid point:
1import numpy 2from pysdic import compute_shape_functions 3 4coords = numpy.array([[0.0, 0.0], [0.5, 0.0], [0.0, 0.5], [0.7, 0.5]]) 5shape_funcs = compute_shape_functions(coords, element_type='triangle_6') 6print("Shape function values for triangle_6:") 7print(shape_funcs)
Shape function values for triangle_6: [[ 1. 0. 0. 0. 0. 0.] [ 0. 0. 0. 1. 0. 0.] [ 0. 0. 0. 0. 1. 0.] [ 0. 0. 0. 0. 0. 0.]]