Skip to content

Commit 00eec7d

Browse files
authored
Merge pull request #289 from rubenwiersma/main
Add Hessian and Hessian energy functions to Python bindings
2 parents 0fbf1c9 + a6cc016 commit 00eec7d

File tree

4 files changed

+144
-0
lines changed

4 files changed

+144
-0
lines changed

src/curved_hessian_energy.cpp

Lines changed: 46 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,46 @@
1+
#include "default_types.h"
2+
#include <igl/curved_hessian_energy.h>
3+
#include <nanobind/nanobind.h>
4+
#include <nanobind/ndarray.h>
5+
#include <nanobind/eigen/dense.h>
6+
#include <nanobind/eigen/sparse.h>
7+
#include <nanobind/stl/tuple.h>
8+
9+
namespace nb = nanobind;
10+
using namespace nb::literals;
11+
12+
namespace pyigl
13+
{
14+
auto curved_hessian_energy(
15+
const nb::DRef<const Eigen::MatrixXN> &V,
16+
const nb::DRef<const Eigen::MatrixXI> &F)
17+
{
18+
Eigen::SparseMatrixN Q;
19+
igl::curved_hessian_energy(V,F,Q);
20+
return Q;
21+
}
22+
23+
}
24+
25+
// Bind the wrapper to the Python module
26+
void bind_curved_hessian_energy(nb::module_ &m)
27+
{
28+
m.def(
29+
"curved_hessian_energy",
30+
&pyigl::curved_hessian_energy,
31+
"V"_a,
32+
"F"_a,
33+
R"(Computes the curved Hessian energy using the Crouzeix-Raviart discretization.
34+
See Oded Stein, Alec Jacobson, Max Wardetzky, Eitan Grinspun, 2020.
35+
"A Smoothness Energy without Boundary Distortion for Curved Surfaces"
36+
37+
@tparam DerivedV derived type of eigen matrix for V (e.g. derived from
38+
MatrixXd)
39+
@tparam DerivedF derived type of eigen matrix for F (e.g. derived from
40+
MatrixXi)
41+
@tparam Scalar scalar type for eigen sparse matrix (e.g. double)
42+
@param[in] V #V by dim list of mesh vertex positions
43+
@param[in] F #F by 3 list of mesh elements (must be triangles)
44+
@param[out] Q #V by #V Hessian energy matrix, each row/column i corresponding to V(i,:)");
45+
46+
}

src/hessian.cpp

Lines changed: 48 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,48 @@
1+
#include "default_types.h"
2+
#include <igl/hessian.h>
3+
#include <nanobind/nanobind.h>
4+
#include <nanobind/ndarray.h>
5+
#include <nanobind/eigen/dense.h>
6+
#include <nanobind/eigen/sparse.h>
7+
#include <nanobind/stl/tuple.h>
8+
9+
namespace nb = nanobind;
10+
using namespace nb::literals;
11+
12+
namespace pyigl
13+
{
14+
auto hessian(
15+
const nb::DRef<const Eigen::MatrixXN> &V,
16+
const nb::DRef<const Eigen::MatrixXI> &F)
17+
{
18+
Eigen::SparseMatrixN H;
19+
igl::hessian(V,F,H);
20+
return H;
21+
}
22+
23+
}
24+
25+
// Bind the wrapper to the Python module
26+
void bind_hessian(nb::module_ &m)
27+
{
28+
m.def(
29+
"hessian",
30+
&pyigl::hessian,
31+
"V"_a,
32+
"F"_a,
33+
R"(Constructs the finite element Hessian matrix
34+
as described in https://arxiv.org/abs/1707.04348,
35+
Natural Boundary Conditions for Smoothing in Geometry Processing
36+
(Oded Stein, Eitan Grinspun, Max Wardetzky, Alec Jacobson)
37+
The interior vertices are NOT set to zero yet.
38+
39+
@tparam DerivedV derived type of eigen matrix for V (e.g. derived from
40+
MatrixXd)
41+
@tparam DerivedF derived type of eigen matrix for F (e.g. derived from
42+
MatrixXi)
43+
@tparam Scalar scalar type for eigen sparse matrix (e.g. double)
44+
@param[in] V #V by dim list of mesh vertex positions
45+
@param[in] F #F by 3 list of mesh elements (must be triangles)
46+
@param[out] H dim²⋅#V by #V Hessian matrix, each column i corresponding to V(i,:)");
47+
48+
}

src/hessian_energy.cpp

Lines changed: 47 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,47 @@
1+
#include "default_types.h"
2+
#include <igl/hessian_energy.h>
3+
#include <nanobind/nanobind.h>
4+
#include <nanobind/ndarray.h>
5+
#include <nanobind/eigen/dense.h>
6+
#include <nanobind/eigen/sparse.h>
7+
#include <nanobind/stl/tuple.h>
8+
9+
namespace nb = nanobind;
10+
using namespace nb::literals;
11+
12+
namespace pyigl
13+
{
14+
auto hessian_energy(
15+
const nb::DRef<const Eigen::MatrixXN> &V,
16+
const nb::DRef<const Eigen::MatrixXI> &F)
17+
{
18+
Eigen::SparseMatrixN Q;
19+
igl::hessian_energy(V,F,Q);
20+
return Q;
21+
}
22+
23+
}
24+
25+
// Bind the wrapper to the Python module
26+
void bind_hessian_energy(nb::module_ &m)
27+
{
28+
m.def(
29+
"hessian_energy",
30+
&pyigl::hessian_energy,
31+
"V"_a,
32+
"F"_a,
33+
R"(Constructs the Hessian energy matrix using mixed FEM
34+
as described in https://arxiv.org/abs/1707.04348
35+
Natural Boundary Conditions for Smoothing in Geometry Processing
36+
(Oded Stein, Eitan Grinspun, Max Wardetzky, Alec Jacobson)
37+
38+
@tparam DerivedV derived type of eigen matrix for V (e.g. derived from
39+
MatrixXd)
40+
@tparam DerivedF derived type of eigen matrix for F (e.g. derived from
41+
MatrixXi)
42+
@tparam Scalar scalar type for eigen sparse matrix (e.g. double)
43+
@param[in] V #V by dim list of mesh vertex positions
44+
@param[in] F #F by 3 list of mesh elements (must be triangles)
45+
@param[out] Q #V by #V Hessian energy matrix, each row/column i corresponding to V(i,:)");
46+
47+
}

tests/test_all.py

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -109,6 +109,9 @@ def test_operators():
109109
K = igl.cotmatrix_entries(l=l)
110110
B1,B2,B3 = igl.local_basis(V,F)
111111
G = igl.grad(V,F)
112+
H = igl.hessian(V, F)
113+
QH = igl.hessian_energy(V, F)
114+
cQH = igl.curved_hessian_energy(V, F)
112115

113116
M = igl.massmatrix(V,F)
114117
M = igl.massmatrix(V,F,type=igl.MASSMATRIX_TYPE_BARYCENTRIC)

0 commit comments

Comments
 (0)