pysdic.Mesh.compute_shape_functions#

Mesh.compute_shape_functions(integration_points, return_derivatives=False, precompute=False, default=0.0)[source]#

Compute the shape functions and their first derivatives for the mesh elements at given IntegrationPoints (convenient method for pysdic.compute_shape_functions()). The element_type must be specified in the mesh connectivity to compute the shape functions for the elements of the mesh.

In a space of dimension \(E\), we consider a \(K\)-dimensional element (with \(K \leq E\)) defined by \(N_{vpe}\) nodes/vertices.

For an \(K\)-dimensional element defined by \(N_{vpe}\) nodes, points inside the element are represented in a local coordinate system \((\xi, \eta, \zeta, ...)\) also named natural coordinates. Shape functions are defined in this local coordinate system in order to interpolate values at any point within the element based on the values at the nodes.

\[P(\xi, \eta, \zeta, ...) = \sum_{i=1}^{N_{vpe}} N_i(\xi, \eta, \zeta, ...) P_i\]

where \(P\) is the interpolated value at the point, \(N_i\) are the shape functions, and \(P_i\) are the nodal values.

Note

The function is a convenience wrapper around the pysdic.compute_shape_functions() function.

Warning

No test are performed to check if the provided integration points are consistent with the mesh (i.e., existing element indices). Only shape check is performed. The behavior of the function is undefined if the provided integration points are not consistent with the mesh.

Parameters:
  • integration_points (IntegrationPoints) – An instance of the IntegrationPoints class containing the natural coordinates of the points where the shape functions are to be evaluated.

  • return_derivatives (bool, optional) – If set to True, the method will also compute and return the derivatives of the shape functions with respect to the local coordinates. Default is False.

  • precompute (bool, optional) – If set to True, the computed shape functions (and their derivatives if return_derivatives is True) will be stored in the integration points object for future reuse, avoiding redundant computations if the shape functions are needed again at the same integration points under the keys shape_functions and shape_function_derivatives. Default is False.

  • default (Real, optional) – The default value to assign to shape functions for points outside the valid range. Default is 0.0.

Returns:

  • shape_functions (numpy.ndarray) – An array of shape (\(M\), \(N_{vpe}\)) containing the evaluated shape functions at the specified points for the \(N_{vpe}\) nodes of the element.

  • shape_function_derivatives (numpy.ndarray, optional) – An array of shape (\(M\), \(N_{vpe}\), \(K\)) containing the derivatives of the shape functions with respect to the local coordinates, if return_derivatives is True.

Return type:

ndarray | Tuple[ndarray, ndarray]

See also

pysdic.compute_shape_functions()

For more information on the algorithm used to compute the shape functions and their derivatives, and an example of how to use it.