pysdic.assemble_jacobian_matrix#

assemble_jacobian_matrix(shape_function_derivatives, remapped_coordinates)[source]#

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)\) from precomputed shape function derivatives and remapped vertex coordinates.

jacobians[i] contains the Jacobian matrix for the integration point \(i\).

Note

  • Inputs shape_function_derivatives and remapped_coordinates will be converted to numpy.float64.

Warning

No tests are performed to check if the inputs array are consistent (e.g., if the connectivity contains invalid vertex indices, …). The behavior of the function is undefined in this case. Only shapes are validated.

Parameters:
  • shape_function_derivatives (ArrayLike) – An array of shape \((N_{p}, N_{vpe}, K)\) containing the derivatives of the shape functions with respect to the local coordinates, evaluated at \(N_{p}\) points for the \(N_{vpe}\) nodes of the element.

  • remapped_coordinates (ArrayLike) – An array of shape \((N_{p}, N_{vpe}, E)\) containing the remapped global coordinates of the vertices for each integration point within the elements

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:

numpy.ndarray

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.

\[\nabla_{\Xi} N = J^{T} \nabla_{X} N = J_{X/\Xi}^{T} \nabla_{X} N\]

To obtain the jacobian matrix for the transformation from global coordinates to natural coordinates noted \(J_{\Xi/X}\) such as:

\[\nabla_{X} N = J_{\Xi/X}^{T} \nabla_{\Xi} N\]

we get :

\[\nabla_{X} N = (J^{T})^{-1} \nabla_{\Xi} N = (J^{-1})^{T} \nabla_{\Xi} N\]

So the jacobian matrix for the transformation from global coordinates to natural coordinates is the inverse of the jacobian matrix for the transformation from natural coordinates to global coordinates.

\[J_{\Xi/X} = J^{-1} = J_{X/\Xi}^{-1}\]

For \(K \le E\), the Jacobian matrix is not square and thus not invertible. In this case, the pseudo-inverse can be used to compute the transformation from global coordinates to natural coordinates. The Moore-Penrose pseudo-inverse of the Jacobian matrix \(J \in \mathbb{R}^{E \times K}\) is given by \(J^{+} = (J^{T} J)^{-1} J^{T}\) because \(J\) has full column rank (i.e., the columns of \(J\) are linearly independent).

See also

pysdic.compute_shape_functions

To compute the shape functions and their derivatives at given natural coordinates within elements.

pysdic.remap_vertices_coordinates

To remap the vertices coordinates for each element at given integration points.

pysdic.compute_jacobian_matrix

To compute the Jacobian matrices directly from natural coordinates and mesh connectivity.

Examples

Lets construct a simple 2D mesh and compute the Jacobian at given integration points for triangular elements.

 1import numpy
 2from pysdic import (
 3    remap_vertices_coordinates,
 4    assemble_jacobian_matrix,
 5    compute_triangle_3_shape_functions
 6)
 7
 8vertices_coordinates = numpy.array(
 9    [[0.0, 0.0],
10     [1.0, 0.0],
11     [1.0, 1.0],
12     [0.0, 1.0]]
13)
14
15connectivity = numpy.array(
16    [[0, 1, 2],
17     [0, 2, 3]]
18)
19
20natural_coordinates = numpy.array(
21    [[0.2, 0.3],
22     [0.6, 0.2]]
23)
24
25element_indices = numpy.array([0, 1])
26
27remapped_coords = remap_vertices_coordinates(
28    vertices_coordinates=vertices_coordinates,
29    connectivity=connectivity,
30    element_indices=element_indices
31)
32
33_, shape_func_derivs = compute_triangle_3_shape_functions(
34    natural_coordinates=natural_coordinates,
35    return_derivatives=True
36)
37
38jacobians = assemble_jacobian_matrix(
39    shape_function_derivatives=shape_func_derivs,
40    remapped_coordinates=remapped_coords
41)
42
43print(f"jacobians (shape={jacobians.shape}):")
44print(jacobians)
jacobians (shape=(2, 2, 2)):
[[[1. 0.]
  [1. 1.]]]

 [[1. 1.]
  [0. 1.]]]