Skip to content

Commit 7914fbe

Browse files
Claus SusanneClaus Susanne
Claus Susanne
authored and
Claus Susanne
committed
add test to use a struct in C-function similar to tabulate_tensor using void*
1 parent abe2391 commit 7914fbe

File tree

1 file changed

+114
-0
lines changed

1 file changed

+114
-0
lines changed

test/test_custom_data.py

+114
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,114 @@
1+
# Copyright (C) 2024 Susanne Claus
2+
#
3+
# This file is part of FFCx. (https://www.fenicsproject.org)
4+
#
5+
# SPDX-License-Identifier: LGPL-3.0-or-later
6+
7+
import numpy as np
8+
from cffi import FFI
9+
10+
# Define custom tabulate tensor function in C with a struct
11+
# Step 1: Define the function in C and set up the CFFI builder
12+
ffibuilder = FFI()
13+
ffibuilder.set_source(
14+
"_cffi_kernelA",
15+
r"""
16+
typedef struct {
17+
uint8_t size;
18+
double* values;
19+
} cell_data;
20+
21+
void tabulate_tensor_integral_add_values(double* restrict A,
22+
const double* restrict w,
23+
const double* restrict c,
24+
const double* restrict coordinate_dofs,
25+
const int* restrict entity_local_index,
26+
const uint8_t* restrict quadrature_permutation,
27+
void* custom_data)
28+
{
29+
// Cast the void* custom_data to cell_data*
30+
cell_data* custom_data_ptr = (cell_data*)custom_data;
31+
32+
// Access the custom data
33+
uint8_t size = custom_data_ptr->size;
34+
double* values = custom_data_ptr->values;
35+
36+
// Use the values in your computations
37+
for (uint8_t i = 0; i < size; i++) {
38+
A[0] += values[i];
39+
}
40+
}
41+
""",
42+
)
43+
ffibuilder.cdef(
44+
"""
45+
typedef struct {
46+
uint8_t size;
47+
double* values;
48+
} cell_data;
49+
50+
void tabulate_tensor_integral_add_values(double* restrict A,
51+
const double* restrict w,
52+
const double* restrict c,
53+
const double* restrict coordinate_dofs,
54+
const int* restrict entity_local_index,
55+
const uint8_t* restrict quadrature_permutation,
56+
void* custom_data);
57+
"""
58+
)
59+
60+
# Step 2: Compile the C code
61+
ffibuilder.compile(verbose=True)
62+
63+
64+
def test_tabulate_tensor_integral_add_values():
65+
# Step 3: Import the compiled library
66+
from _cffi_kernelA import ffi, lib
67+
68+
# Define cell data
69+
size = 2
70+
values = np.array([2.0, 1.0], dtype=np.float64)
71+
expected_result = np.array([3.0], dtype=np.float64)
72+
73+
# Define the input arguments
74+
A = np.zeros(1, dtype=np.float64)
75+
w = np.array([1.0], dtype=np.float64)
76+
c = np.array([0.0], dtype=np.float64)
77+
coordinate_dofs = np.array(
78+
[0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, 0.0], dtype=np.float64
79+
)
80+
entity_local_index = np.array([0], dtype=np.int32)
81+
quadrature_permutation = np.array([0], dtype=np.uint8)
82+
83+
# Cast the arguments to the appropriate C types
84+
A_ptr = ffi.cast("double*", A.ctypes.data)
85+
w_ptr = ffi.cast("double*", w.ctypes.data)
86+
c_ptr = ffi.cast("double*", c.ctypes.data)
87+
coordinate_dofs_ptr = ffi.cast("double*", coordinate_dofs.ctypes.data)
88+
entity_local_index_ptr = ffi.cast("int*", entity_local_index.ctypes.data)
89+
quadrature_permutation_ptr = ffi.cast("uint8_t*", quadrature_permutation.ctypes.data)
90+
91+
# Use ffi.from_buffer to create a CFFI pointer from the NumPy array
92+
values_ptr = ffi.from_buffer(values)
93+
94+
# Allocate memory for the struct
95+
custom_data = ffi.new("cell_data*")
96+
custom_data.size = size
97+
custom_data.values = values_ptr
98+
99+
# Cast the struct to void*
100+
custom_data_ptr = ffi.cast("void*", custom_data)
101+
102+
# Call the function
103+
lib.tabulate_tensor_integral_add_values(
104+
A_ptr,
105+
w_ptr,
106+
c_ptr,
107+
coordinate_dofs_ptr,
108+
entity_local_index_ptr,
109+
quadrature_permutation_ptr,
110+
custom_data_ptr,
111+
)
112+
113+
# Assert the result
114+
np.testing.assert_allclose(A, expected_result, rtol=1e-5)

0 commit comments

Comments
 (0)