Skip to content

Commit

Permalink
#24 Added first (faulty) implementation of 1-bit per bool morphology.
Browse files Browse the repository at this point in the history
  • Loading branch information
carljohnsen committed May 14, 2024
1 parent 43bdb3d commit fd9c463
Show file tree
Hide file tree
Showing 3 changed files with 149 additions and 11 deletions.
102 changes: 102 additions & 0 deletions src/lib/cpp/cpu_seq/morphology.cc
Original file line number Diff line number Diff line change
Expand Up @@ -46,4 +46,106 @@ void morphology_3d_sphere(
}
}

template <typename Op, uint32_t neutral>
void morphology_3d_sphere_bitpacked(
const uint32_t *voxels,
const int64_t radius,
const int64_t N[3],
const int64_t strides[3],
uint32_t *result) {
// TODO assumes that Nx is a multiple of 32, which is true for scale <= 4
Op op;
int64_t
k = radius*2 + 1,
sqradius = radius * radius;

// TODO handle k < 32
// TODO templated construction? Has to 'hardcode' a radius, but that is beneficial anyways.
// Create the kernel
uint32_t *kernel = (uint32_t*) malloc(k*k*sizeof(uint32_t));

#pragma omp parallel for collapse(2)
for (int64_t z = -radius; z <= radius; z++) {
for (int64_t y = -radius; y <= radius; y++) {
uint32_t row = 0;
for (int64_t x = 0; x < 32; x++) {
uint32_t element = (x-radius)*(x-radius) + y*y + z*z <= sqradius;
row |= element << (31 - x);
}
kernel[(z+radius)*k + y+radius] = row;
//printf("%08x ", row);
}
//printf("\n");
}

#pragma omp parallel for collapse(3)
for (int64_t z = 0; z < N[0]; z++) {
for (int64_t y = 0; y < N[1]; y++) {
for (int64_t x = 0; x < N[2]; x++) {
// Compute boundaries
int64_t flat_index = z*strides[0] + y*strides[1] + x*strides[2];
int64_t X[3] = {z, y, x};
int64_t limits[6];
for (int axis = 0; axis < 3; axis++) {
limits[(axis*2)] = -min(radius, X[axis]);
limits[(axis*2)+1] = min(radius, N[axis] - X[axis] - 1);
}

// Apply the spherical kernel
uint32_t value = neutral;
for (int64_t pz = limits[0]; pz <= limits[1]; pz++) {
for (int64_t py = limits[2]; py <= limits[3]; py++) {
uint32_t
voxels_row = voxels[flat_index / 32],
kernel_row = kernel[pz*k + py];

int64_t
beginning = x - radius,
end = x + radius;
if (beginning < 0) { // Case 1
int64_t
mask_shift = std::abs(beginning),
neutral_shift = 32 - mask_shift;
//if (z == 16 && y == 16) printf("Case 1: %ld %08x %08x %08x %08x\n", x, voxels_row, kernel_row, voxels_row >> mask_shift, neutral << neutral_shift);
voxels_row = voxels_row >> mask_shift | neutral << neutral_shift;
} else if (beginning / 32 != x / 32) { // Case 2
int64_t
mask1_shift = beginning % 32,
mask0_shift = 32 - mask1_shift;
uint32_t voxels1 = voxels[(flat_index / 32) - 1];
voxels_row = (voxels1 << mask1_shift) | (voxels_row >> mask0_shift);
} else if (beginning / 32 == end / 32) { // Case 3
int64_t
mask_shift = beginning % 32,
neutral_shift = 32 - mask_shift;
voxels_row = voxels_row << mask_shift | neutral >> neutral_shift;
} else if(end / 32 != x / 32) { // Case 4
int64_t
mask0_shift = beginning % 32,
mask1_shift = 32 - mask0_shift;
uint32_t voxels1 = voxels[(flat_index / 32) + 1];
voxels_row = (voxels_row << mask0_shift) | (voxels1 >> mask1_shift);
} else if (end >= N[2]) { // Case 5
int64_t
mask_shift = beginning % 32,
neutral_shift = 32 - mask_shift;
voxels_row = voxels_row << mask_shift | neutral >> neutral_shift;
} else {
assert (false && "Should not reach this point - some case is missing.");
}

value = op(value, voxels_row & kernel_row);
}
}
// dilate:
value = (value != 0) << (31 - x % 32);
// erode:

// Store the results
result[flat_index/32] |= value;
}
}
}
}

} // namespace cpu_seq
8 changes: 8 additions & 0 deletions src/lib/cpp/include/morphology.hh
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,14 @@ void morphology_3d_sphere(
const int64_t strides[3],
mask_type *result);

template <typename Op, uint32_t neutral>
void morphology_3d_sphere_bitpacked(
const uint32_t *voxels,
const int64_t radius,
const int64_t N[3],
const int64_t strides[3],
uint32_t *result);

} // namespace NS

#endif
50 changes: 39 additions & 11 deletions src/pybind/morphology-pybind.cc
Original file line number Diff line number Diff line change
Expand Up @@ -26,37 +26,65 @@ void morphology_3d_sphere_wrapper(
mask_type *result = static_cast<mask_type*>(result_info.ptr);

if (radius == (int64_t) 16) {
NS::morphology_3d_sphere_r16<Op, neutral>(voxels, N, strides, result);
//NS::morphology_3d_sphere_r16<Op, neutral>(voxels, N, strides, result);
} else {
NS::morphology_3d_sphere<Op, neutral>(voxels, radius, N, strides, result);
}
}

template <typename Op, bool neutral>
void morphology_3d_sphere_wrapper_alt(
const py::array_t<mask_type> &np_voxels,
template <typename Op, uint32_t neutral>
void morphology_3d_sphere_bitpacked(
const uint32_t *voxels,
const int64_t radius,
py::array_t<mask_type> np_result) {
const int64_t N[3],
const int64_t strides[3],
uint32_t *result);

template <typename Op, uint32_t neutral>
void morphology_3d_sphere_bitpacked_wrapper(
const py::array_t<uint32_t> &np_voxels,
const int64_t radius,
py::array_t<uint32_t> np_result) {
auto
voxels_info = np_voxels.request(),
result_info = np_result.request();

int64_t Nz = voxels_info.shape[0], Ny = voxels_info.shape[1], Nx = voxels_info.shape[2];
int64_t Nz = voxels_info.shape[0], Ny = voxels_info.shape[1], Nx = voxels_info.shape[2] * 32;
int64_t N[3] = {Nz, Ny, Nx};
int64_t strides[3] = {Ny*Nx, Nx, 1};

const mask_type *voxels = static_cast<const mask_type*>(voxels_info.ptr);
mask_type *result = static_cast<mask_type*>(result_info.ptr);
const uint32_t *voxels = static_cast<const uint32_t*>(voxels_info.ptr);
uint32_t *result = static_cast<uint32_t*>(result_info.ptr);

NS::morphology_3d_sphere_alt(voxels, radius, N, strides, result);
NS::morphology_3d_sphere_bitpacked<Op, neutral>(voxels, radius, N, strides, result);
}

//template <typename Op, bool neutral>
//void morphology_3d_sphere_wrapper_alt(
// const py::array_t<mask_type> &np_voxels,
// const int64_t radius,
// py::array_t<mask_type> np_result) {
// auto
// voxels_info = np_voxels.request(),
// result_info = np_result.request();
//
// int64_t Nz = voxels_info.shape[0], Ny = voxels_info.shape[1], Nx = voxels_info.shape[2];
// int64_t N[3] = {Nz, Ny, Nx};
// int64_t strides[3] = {Ny*Nx, Nx, 1};
//
// const mask_type *voxels = static_cast<const mask_type*>(voxels_info.ptr);
// mask_type *result = static_cast<mask_type*>(result_info.ptr);
//
// NS::morphology_3d_sphere_alt(voxels, radius, N, strides, result);
//}

} // namespace python_api

PYBIND11_MODULE(morphology, m) {
m.doc() = "Morphology operations."; // optional module docstring
m.def("dilate_3d_sphere", &python_api::morphology_3d_sphere_wrapper<std::bit_or<mask_type>, false>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
m.def("dilate_3d_sphere_alt", &python_api::morphology_3d_sphere_wrapper_alt<std::bit_or<mask_type>, false>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
m.def("dilate_3d_sphere_bitpacked", &python_api::morphology_3d_sphere_bitpacked_wrapper<std::bit_or<uint32_t>, 0>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
//m.def("dilate_3d_sphere_alt", &python_api::morphology_3d_sphere_wrapper_alt<std::bit_or<mask_type>, false>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
m.def("erode_3d_sphere", &python_api::morphology_3d_sphere_wrapper<std::bit_and<mask_type>, true>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
m.def("erode_3d_sphere_alt", &python_api::morphology_3d_sphere_wrapper_alt<std::bit_and<mask_type>, true>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
//m.def("erode_3d_sphere_alt", &python_api::morphology_3d_sphere_wrapper_alt<std::bit_and<mask_type>, true>, py::arg("np_voxels"), py::arg("radius"), py::arg("np_result").noconvert());
}

0 comments on commit fd9c463

Please sign in to comment.