pysdic.assemble_property_interpolation#
- assemble_property_interpolation(property_array, connectivity, element_indices, shape_functions, *, skip_m1=True, default=nan)[source]#
Interpolate a property defined at the nodes of a mesh to given integration points within elements using precomputed shape functions .
interpolated_properties[i]contains the interpolated property at the integration point \(i\).Note
Inputs
property_arrayandshape_functionswill 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:
property_array (ArrayLike) – An array of shape \((N_{v}, P)\) containing the property values defined at the nodes of the mesh. If 1D-array is provided, it will be treated as a single-component property of shape \((N_{v}, 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_indices (ArrayLike) – An array of shape \((N_{p},)\) containing the indices of each element corresponding to the \(N_{p}\) integration points.
shape_functions (ArrayLike) – An array of shape \((N_{p}, N_{vpe})\) containing the shape function values evaluated at \(N_{p}\) points for the \(N_{vpe}\) nodes of the element.
skip_m1 (
bool, optional) – If set toTrue, any element index of -1 inelement_indiceswill result in the corresponding interpolated property being set todefault. Default isTrue.default (Union[Real, ArrayLike], optional) – The default value to assign to interpolated properties for integration points associated with an element index of -1 when
skip_m1isTrue. Default isnumpy.nan. The input can also be anumpy.ndarrayof shape (P,) to assign different default values for each property component.
- Returns:
interpolated_properties – An array of shape \((N_{p}, P)\) containing the interpolated property values at each of the \(N_{p}\) integration points.
- 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.
For a given set of \(N_{p}\) integration points located within elements, the interpolated property array has shape \((N_{p}, P)\) where \(P\) is the number of property components (e.g., 1 for scalar properties, 3 for vector properties).
The property at each integration point is computed as:
\[P(\xi, \eta, \zeta, ...) = \sum_{i=1}^{N_{vpe}} N_i(\xi, \eta, \zeta, ...) P_i\]where \(N_i\) are the shape functions associated with each node, and \(P_i\) are the property values at the nodes of the element.
See also
pysdic.compute_shape_functionsTo compute the shape functions at given natural coordinates within elements.
pysdic.compute_property_interpolationTo compute the property interpolation directly from mesh and integration point information.
pysdic.assemble_property_interpolationTo project properties defined at integration points back to the nodes of the mesh.
Examples
Lets construct a simple 2D mesh and interpolate a scalar property at given integration points for triangular elements.
1import numpy 2from pysdic import ( 3 compute_shape_functions, 4 assemble_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 19natural_coordinates = numpy.array( 20 [[0.2, 0.3], 21 [0.6, 0.2]] 22) 23 24element_indices = numpy.array([0, 1]) 25 26property_array = numpy.array([10.0, 20.0, 30.0, 40.0]) # Scalar property 27 28shape_functions = compute_shape_functions( 29 natural_coordinates=natural_coordinates, 30 element_type='triangle_3', 31 return_derivatives=False 32) 33 34interpolated_props = assemble_property_interpolation( 35 property_array=property_array, 36 connectivity=connectivity, 37 element_indices=element_indices, 38 shape_functions=shape_functions 39) 40 41print(f"interpolated properties (shape={interpolated_props.shape}):") 42print(interpolated_props)
interpolated properties (shape=(2, 1)): [[18. ] [28. ]]