pysdic.compute_jacobian_matrix#
- compute_jacobian_matrix(vertices_coordinates, connectivity, element_type, natural_coordinates, element_indices, *, skip_m1=True, default=nan)[source]#
Compute and assemble the Jacobian (derivative) matrix \(J\) of the transformation from natural coordinates \((\xi, \eta, \zeta, ...)\) to global coordinates \((x, y, z, ...)\) with shape \((N_{p}, E, K)\) at given integration points within specified elements.
jacobians[i]contains the Jacobian matrix for the integration point \(i\).This function combines the remapping of vertex coordinates, computation of shape function derivatives, and assembly of the Jacobian matrices for the entire mesh.
Note
Inputs
vertices_coordinatesandnatural_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:
vertices_coordinates (ArrayLike) – An array of shape \((N_{v}, E)\) containing the global coordinates of the vertices in the mesh.
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.
skip_m1 (
bool, optional) – If set toTrue, any element index of -1 inelement_indiceswill result in the corresponding Jacobian being set todefault. Default isTrue.default (Real, optional) – The default value to assign to Jacobians for integration points associated with an element index of -1 when
skip_m1isTrue. The out-of-range natural coordinates will also result in Jacobians being set to this value. Default isnumpy.nan.
- Returns:
jacobians – An array of shape \((N_{p}, E, K)\) containing the Jacobian matrices for each of the \(N_{p}\) points. Each Jacobian matrix has dimensions \((E \times K)\).
- Return type:
Notes
In a space of dimension \(E\) with a mesh constituted of \(N_{e}\) elements and \(N_{v}\) nodes, the mesh is composed of \(K\)-dimensional elements (with \(K \leq E\)) defined by \(N_{vpe}\) nodes for each element.
The Jacobian matrix has dimensions \((E \times K)\) and constructed as follows:
\[\begin{split}J = J_{X/\Xi} = \frac{dX}{d\Xi} = \begin{bmatrix} \frac{\partial x}{\partial \xi} & \frac{\partial x}{\partial \eta} & \frac{\partial x}{\partial \zeta} & \ldots \\ \frac{\partial y}{\partial \xi} & \frac{\partial y}{\partial \eta} & \frac{\partial y}{\partial \zeta} & \ldots \\ \frac{\partial z}{\partial \xi} & \frac{\partial z}{\partial \eta} & \frac{\partial z}{\partial \zeta} & \ldots \\ \vdots & \vdots & \vdots & \ddots \end{bmatrix}\end{split}\]where each entry is computed as:
\[\frac{\partial x_i}{\partial \xi_j} = \sum_{a=1}^{N_{vpe}} \frac{\partial N_a}{\partial \xi_j} x_{i,a}\]Here, \(N_a\) are the shape functions associated with each node, and \(x_{i,a}\) are the global coordinates of the nodes.
See the notes in
pysdic.assemble_jacobian_matrix()for more details on the transformation between natural and global coordinates and the use of pseudo-inverse for non-square Jacobian matrices.See also
pysdic.compute_shape_functionsTo compute the shape functions and their derivatives at given natural coordinates within elements.
pysdic.remap_vertices_coordinatesTo remap the vertices coordinates for each element at given integration points.
pysdic.assemble_jacobian_matrixTo assemble the Jacobian matrices from shape function derivatives and remapped coordinates.
Examples
Lets construct a simple 2D mesh and compute the Jacobian at given integration points for triangular elements.
1import numpy 2from pysdic import compute_jacobian_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 23jacobians = compute_jacobian_matrix( 24 vertices_coordinates=vertices_coordinates, 25 connectivity=connectivity, 26 element_type='triangle_3', 27 natural_coordinates=natural_coordinates, 28 element_indices=element_indices, 29) 30 31print(f"jacobians (shape={jacobians.shape}):") 32print(jacobians)
jacobians (shape=(2, 2, 2)): [[[1. 1.] [0. 1.]] [[1. 0.] [1. 1.]]]