pysdic.compute_property_projection#

compute_property_projection(property_array, connectivity, element_type, natural_coordinates, element_indices, point_weights=None, n_vertices=None, *, sparse=False, skip_m1=True, return_unaffected=False)[source]#

Project a property defined at integration points within elements back to the nodes of a mesh using mesh and integration point information.

This function combines the remapping of vertex coordinates, computation of shape functions, and assembly of the property projection for the entire mesh.

projected_properties[i] contains the projected property at the node \(i\).

Note

  • Inputs property_array, natural_coordinates and point_weights will be converted to numpy.float64.

  • Inputs connectivity and element_indices will be converted to numpy.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 -1 in element_indices for invalid elements, ensure to set skip_m1 to True to avoid indexing errors.

Parameters:
  • property_array (ArrayLike) – An array of shape \((N_{p}, P)\) containing the property values defined at the \(N_{p}\) integration points. If 1D-array is provided, it will be treated as a single-component property of shape \((N_{p}, 1)\).

  • 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.

  • point_weights (Optional[ArrayLike], optional) – An array of shape \((N_{p},)\) containing the weights associated with each integration point. If not provided, all weights will be assumed to be equal to one.

  • 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 index in connectivity plus one.

  • sparse (bool, optional) – If set to True, the shape functions matrix will be constructed as a sparse matrix to optimize memory usage for large meshes. Default is False.

  • skip_m1 (bool, optional) – If set to True, any element index of -1 in element_indices will result in the corresponding integration point being ignored during the projection. Default is True.

  • return_unaffected (bool, optional) – If set to True, the function will return a tuple containing the projected properties and a boolean mask indicating which nodes were unaffected by any integration point. Default is False.

Returns:

  • projected_properties (numpy.ndarray) – An array of shape \((N_{v}, P)\) containing the projected property values at the nodes of the mesh.

  • unaffected_mask (numpy.ndarray, optional) – A boolean array of shape \((N_{v},)\) indicating which nodes were unaffected by any integration point. This is only returned if return_unaffected is set to True. True indicates the node was unaffected, while False indicates it was affected. So using projected_properties[unaffected_mask, :] = default will set the unaffected nodes to default.

Return type:

ndarray | Tuple[ndarray, ndarray]

Notes

See the documentation of assemble_property_projection() for the mathematical formulation and demonstration of the projection process.

See also

pysdic.compute_shape_functions

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

pysdic.assemble_property_projection

To project properties defined at integration points back to the nodes of the mesh using a precomputed shape functions matrix.

pysdic.compute_property_interpolation

To interpolate properties defined at the nodes of the mesh to integration points within elements.

Examples

Lets construct a simple 2D mesh and project a scalar property defined at given integration points back to the nodes for triangular elements.

 1import numpy
 2from pysdic import (
 3    compute_property_projection,
 4    compute_property_interpolation
 5)
 6
 7vertices_coordinates = numpy.array(
 8    [[0.0, 0.0],
 9     [1.0, 0.0],
10     [1.0, 1.0],
11     [0.0, 1.0]]
12)
13
14connectivity = numpy.array(
15    [[0, 1, 2],
16     [0, 2, 3]]
17)
18
19N_e = connectivity.shape[0]
20
21property_array = numpy.array([10.0, 20.0, 30.0,  40.0])  # Scalar property
22
23# Define natural coordinates for interpolation (e.g., 4 points per element)
24natural_coordinates = numpy.array(
25    [[0.3, 0.3],
26     [0.2, 0.5],
27     [0.5, 0.2],
28     [0.1, 0.1]]
29)
30N_p = natural_coordinates.shape[0]
31natural_coordinates = numpy.vstack([natural_coordinates] * N_e)
32element_indices = numpy.repeat(numpy.arange(N_e), N_p)
33N_p = natural_coordinates.shape[0]
34
35interpolate_property = compute_property_interpolation(
36    property_array=property_array,
37    connectivity=connectivity,
38    element_type='triangle_3',
39    natural_coordinates=natural_coordinates,
40    element_indices=element_indices
41)
42
43# Now project back to the nodes
44projected_props = compute_property_projection(
45    property_array=interpolate_property,
46    connectivity=connectivity,
47    element_type='triangle_3',
48    natural_coordinates=natural_coordinates,
49    element_indices=element_indices,
50    sparse=False
51)
52
53print(f"projected properties (shape={projected_props.shape}):")
54print(projected_props)
projected properties (shape=(4, 1)):
[[10.]
 [20.]
 [30.]
 [40.]]