Skip to content

Commit 0e2748d

Browse files
committed
fix aabb template nonsense
1 parent 3dcbcbc commit 0e2748d

2 files changed

Lines changed: 31 additions & 18 deletions

File tree

classes/AABB.cpp

Lines changed: 25 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -17,15 +17,26 @@ void check_c_contiguity(py::array_t<T> V, const std::string name = "arg")
1717
throw std::runtime_error(name+" must be C-contiguous with a 2D shape.");
1818
}
1919
};
20+
template <typename T>
21+
void check_f_contiguity(py::array_t<T> V, const std::string name = "arg")
22+
{
23+
auto V_buf = V.request();
24+
if (V_buf.ndim != 2 || V_buf.strides[1] == sizeof(T)) {
25+
throw std::runtime_error(name+" must be F-contiguous (ColMajor) with a 2D shape.");
26+
}
27+
};
2028

2129
template <int DIM>
2230
void init_AABB(py::module_ &m)
2331
{
32+
using namespace Eigen;
2433
using MatrixXdRowMajor = Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic,Eigen::RowMajor>;
2534
using MatrixXiRowMajor = Eigen::Matrix<int,Eigen::Dynamic,Eigen::Dynamic,Eigen::RowMajor>;
26-
using MapMatrixXdRowMajor = Eigen::Map<MatrixXdRowMajor>;
27-
using MapMatrixXiRowMajor = Eigen::Map<MatrixXiRowMajor>;
28-
using DerivedV = MapMatrixXdRowMajor;
35+
using MapMatrixXdRowMajor = Eigen::Map<const MatrixXdRowMajor>;
36+
using MapMatrixXiRowMajor = Eigen::Map<const MatrixXiRowMajor>;
37+
using MapMatrixXd = Eigen::Map<const Eigen::MatrixXd>;
38+
using MapMatrixXi = Eigen::Map<const Eigen::MatrixXi>;
39+
using DerivedV = MapMatrixXd;
2940
using AABB_f64_DIM = igl::AABB<DerivedV,DIM>;
3041
py::class_<AABB_f64_DIM >(m, (std::string("AABB_f64_")+std::to_string(DIM)).c_str() )
3142
.def(py::init([]()
@@ -42,11 +53,11 @@ void init_AABB(py::module_ &m)
4253
throw std::runtime_error("Input matrices must be 2-dimensional.");
4354
}
4455

45-
check_c_contiguity(V,"V");
46-
check_c_contiguity(F,"F");
56+
check_f_contiguity(V,"V");
57+
check_f_contiguity(F,"F");
4758

48-
MapMatrixXdRowMajor V_eigen(static_cast<double*>(V_buf.ptr), V_buf.shape[0], V_buf.shape[1]);
49-
MapMatrixXiRowMajor F_eigen(static_cast<int*>(F_buf.ptr), F_buf.shape[0], F_buf.shape[1]);
59+
MapMatrixXd V_eigen(static_cast<double*>(V_buf.ptr), V_buf.shape[0], V_buf.shape[1]);
60+
MapMatrixXi F_eigen(static_cast<int*>(F_buf.ptr), F_buf.shape[0], F_buf.shape[1]);
5061

5162
self.init(V_eigen, F_eigen);
5263
}, py::arg("V"), py::arg("F"))
@@ -59,20 +70,20 @@ void init_AABB(py::module_ &m)
5970
const bool return_closest_point = false) ->
6071
std::variant<py::object,std::list<py::object> >
6172
{
62-
check_c_contiguity(V,"V");
63-
check_c_contiguity(F,"F");
64-
check_c_contiguity(P,"P");
73+
check_f_contiguity(V,"V");
74+
check_f_contiguity(F,"F");
75+
check_f_contiguity(P,"P");
6576

6677
auto V_buf = V.request();
6778
auto F_buf = F.request();
6879
auto P_buf = P.request();
69-
MapMatrixXdRowMajor V_eigen(static_cast<double*>(V_buf.ptr), V_buf.shape[0], V_buf.shape[1]);
70-
MapMatrixXiRowMajor F_eigen(static_cast<int*>(F_buf.ptr), F_buf.shape[0], F_buf.shape[1]);
71-
MapMatrixXdRowMajor P_eigen(static_cast<double*>(P_buf.ptr), P_buf.shape[0], P_buf.shape[1]);
80+
MapMatrixXd V_eigen(static_cast<double*>(V_buf.ptr), V_buf.shape[0], V_buf.shape[1]);
81+
MapMatrixXi F_eigen(static_cast<int*>(F_buf.ptr), F_buf.shape[0], F_buf.shape[1]);
82+
MapMatrixXd P_eigen(static_cast<double*>(P_buf.ptr), P_buf.shape[0], P_buf.shape[1]);
7283

7384
Eigen::VectorXd sqrD;
7485
Eigen::VectorXi I;
75-
MatrixXdRowMajor C;
86+
MatrixXd C;
7687
self.squared_distance(V_eigen,F_eigen,P_eigen,sqrD,I,C);
7788
if(return_index && return_closest_point)
7889
{

src/in_element.cpp

Lines changed: 6 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -11,7 +11,8 @@
1111
#include <igl/in_element.h>
1212
#include <igl/AABB.h>
1313

14-
using AABB_f64_3 = igl::AABB<Eigen::MatrixXd,3>;
14+
using AABB_f64_3 = igl::AABB<Eigen::Map<const Eigen::MatrixXd>,3>;
15+
using AABB_f64_2 = igl::AABB<Eigen::Map<const Eigen::MatrixXd>,2>;
1516
const char* ds_in_element = R"igl_Qu8mg5v7(
1617
Determine whether each point in a list of points is in the elements of a mesh.
1718
@@ -36,13 +37,13 @@ npe_arg(aabb, AABB_f64_3)
3637
npe_begin_code()
3738

3839
Eigen::VectorXi I;
39-
igl::in_element(V,Ele,Q,aabb,I);
40+
Eigen::Map<const Eigen::MatrixXd> V_map(V.data(),V.rows(),V.cols());
41+
igl::in_element(V_map,Ele,Q,aabb,I);
4042
return npe::move(I);
4143

4244
npe_end_code()
4345

4446

45-
using AABB_f64_2 = igl::AABB<Eigen::MatrixXd,2>;
4647
npe_function(in_element_2)
4748
npe_doc( ds_in_element)
4849
npe_arg(V, Eigen::MatrixXd)
@@ -52,7 +53,8 @@ npe_arg(aabb, AABB_f64_2)
5253
npe_begin_code()
5354

5455
Eigen::VectorXi I;
55-
igl::in_element(V,Ele,Q,aabb,I);
56+
Eigen::Map<const Eigen::MatrixXd> V_map(V.data(),V.rows(),V.cols());
57+
igl::in_element(V_map,Ele,Q,aabb,I);
5658
return npe::move(I);
5759

5860
npe_end_code()

0 commit comments

Comments
 (0)