pysdic.compute_shape_function_matrix#
- compute_shape_function_matrix(connectivity, element_type, natural_coordinates, element_indices, n_vertices=None, *, sparse=False, skip_m1=True, default=0.0)[source]#
Compute and assemble the shape function matrix \(\mathcal{N}\) at given integration points within specified elements with shape \((N_{p}, N_{v})\).
This function combines the computation of shape functions at given natural coordinates with the assembly of the shape function matrix for the entire mesh.
shape_function_matrix[i, j]contains the shape function value at integration point \(i\) for vertex \(j\) and a non-zero value only if vertex \(j\) is part of the element containing integration point \(i\).Note
Input
natural_coordinateswill be converted tonumpy.float64.Inputs
connectivityandelement_indiceswill be converted tonumpy.int64.
Warning
No tests are performed to check if the inputs array are consistent (e.g., if the connectivity contains invalid vertex indices, …). Only shapes are validated. The behavior of the function is undefined in this case.
Important
When using
-1inelement_indicesfor invalid elements, ensure to setskip_m1toTrueto avoid indexing errors.- Parameters:
connectivity (ArrayLike) – An array of shape \((N_{e}, N_{vpe})\) defining the connectivity of the elements in the mesh, where each row contains the indices of the nodes that form an element.
element_type (
str) – A string specifying the type of element (e.g., ‘segment_2’, ‘triangle_3’, etc.) to determine which shape function to use.natural_coordinates (ArrayLike) – An array of shape \((N_{p}, K)\) containing the natural coordinates of the integration points within elements where \(K\) is the topological dimension of the element (e.g., 1 for segments, 2 for triangles/quadrangles, etc.).
element_indices (ArrayLike) – An array of shape \((N_{p},)\) containing the indices of each element corresponding to the \(N_{p}\) integration points.
n_vertices (Optional[Integral], optional) – The total number of vertices \(N_{v}\) in the mesh. If not provided, it will be inferred as the maximum node index in
connectivityplus one.sparse (
bool, optional) – If set toTrue, the function will usescipy.sparseto create a sparse matrix representation of the shape function matrix. Default isFalse.skip_m1 (
bool, optional) – If set toTrue, any element index of -1 inelement_indiceswill result in the corresponding shape function values being set todefault. Default isTrue.default (Real, optional) – The default value to assign to shape function values for integration points associated with an element index of -1 when
skip_m1isTrue. Default is0.0. Only used for non-sparse matrix construction. For sparse matrices, zero-filling is used.
- Returns:
shape_function_matrix – An array of shape \((N_{p}, N_{v})\) or sparse matrix containing the shape function values at each of the \(N_{p}\) integration points for all \(N_{v}\) nodes in the mesh.
- Return type:
Union[numpy.ndarray, scipy.sparse.csr_matrix]
Notes
This function first computes the shape function values at the provided natural coordinates using the appropriate shape function based on the specified
element_type. It then assembles these values into a global shape function matrix for the entire mesh.See also
pysdic.compute_shape_functionsTo compute the shape functions at given natural coordinates within elements with shape \((N_{p}, N_{vpe})\).
pysdic.assemble_shape_function_matrixTo assemble the shape function matrix from precomputed shape function values.
Examples
Lets construct a simple 2D mesh and build the shape function matrix at given integration points for triangular elements.
1import numpy 2from pysdic import compute_shape_function_matrix 3 4vertices_coordinates = numpy.array( 5 [[0.0, 0.0], 6 [1.0, 0.0], 7 [1.0, 1.0], 8 [0.0, 1.0]] 9) 10 11connectivity = numpy.array( 12 [[0, 1, 2], 13 [0, 2, 3]] 14) 15 16natural_coordinates = numpy.array( 17 [[0.2, 0.3], 18 [0.6, 0.2]] 19) 20 21element_indices = numpy.array([0, 1]) 22 23shape_func_matrix = compute_shape_function_matrix( 24 connectivity=connectivity, 25 element_type='triangle_3', 26 natural_coordinates=natural_coordinates, 27 element_indices=element_indices 28) 29 30print(f"shape function matrix (shape={shape_func_matrix.shape}):") 31print(shape_func_matrix)
shape function matrix (shape=(2, 4)): [[0.5 0.2 0.3 0. ] [0.2 0. 0.6 0.2]]