pysdic.compute_segment_2_shape_functions#
- compute_segment_2_shape_functions(natural_coordinates, return_derivatives=False, *, default=0.0)[source]#
Compute the shape functions for a 2-node segment for given
natural_coordinates\(\xi\).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},)\) or \((N_{p}, 1)\), where \(N_{p}\) is the number of points to evaluate.
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 \([-1, 1]\)). By default,
0.0.
- Returns:
shape_functions (
numpy.ndarray) – Shape functions evaluated at the given natural coordinates. The returned array has shape \((N_{p}, 2)\), where each row corresponds to a point and each column to a node.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}, 2, 1)\).
- Return type:
Notes
A 2-node segment represented in the figure below has the following shape functions:
Node No.
\((\xi)\)
Shape Function \(N\)
First Derivative \((\frac{dN}{d\xi})\)
1
\(-1\)
\(N_1(\xi) = \frac{1}{2}(1 - \xi)\)
\(-\frac{1}{2}\)
2
\(1\)
\(N_2(\xi) = \frac{1}{2}(1 + \xi)\)
\(\frac{1}{2}\)
See also
pysdic.compute_segment_3_shape_functionsShape functions for 3-node segment 1D-elements.
pysdic.get_segment_2_gauss_pointsGauss integration points for 2-node segment 1D-elements.
Examples
Compute shape functions without derivatives for 3 valid points and 1 invalid point:
1import numpy 2from pysdic import compute_segment_2_shape_functions 3 4coords = numpy.array([[-1.0], [0.0], [1.0], [1.5]]) 5shape_functions = compute_segment_2_shape_functions(coords) 6print("Shape function values:") 7print(shape_functions)
Shape function values: [[1. 0.] [0.5 0.5] [0. 1.] [0. 0.]]