Skip to content

Commit 7119cb2

Browse files
authored
Merge pull request #273 from sean-d1/main
Add binding for harmonic_integrated_from_laplacian_and_mass
2 parents 1336354 + 1e5af93 commit 7119cb2

File tree

2 files changed

+31
-1
lines changed

2 files changed

+31
-1
lines changed

src/harmonic.cpp

Lines changed: 24 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2,6 +2,7 @@
22
#include <igl/harmonic.h>
33
#include <nanobind/nanobind.h>
44
#include <nanobind/eigen/dense.h>
5+
#include <nanobind/eigen/sparse.h>
56

67
namespace nb = nanobind;
78
using namespace nb::literals;
@@ -23,6 +24,15 @@ namespace pyigl
2324
}
2425
return W;
2526
}
27+
auto harmonic_integrated_from_laplacian_and_mass(
28+
const Eigen::SparseMatrixN &L,
29+
const Eigen::SparseMatrixN &M,
30+
const int k)
31+
{
32+
Eigen::SparseMatrixN Q;
33+
igl::harmonic(L, M, k, Q);
34+
return Q;
35+
}
2636
}
2737

2838
// Bind the wrapper to the Python module
@@ -44,4 +54,17 @@ R"(Compute k-harmonic weight functions "coordinates".
4454
@param[in] bc #b by #W list of boundary values
4555
@param[in] k power of harmonic operation (1: harmonic, 2: biharmonic, etc)
4656
@return W #V by #W list of weights)");
47-
}
57+
m.def(
58+
"harmonic_integrated_from_laplacian_and_mass",
59+
&pyigl::harmonic_integrated_from_laplacian_and_mass,
60+
"L"_a,
61+
"M"_a,
62+
"k"_a,
63+
R"(Build the discrete k-harmonic operator (computing integrated quantities).
64+
That is, if the k-harmonic PDE is Q x = 0, then this minimizes x' Q x.
65+
66+
@param[in] L #V by #V discrete (integrated) Laplacian
67+
@param[in] M #V by #V mass matrix
68+
@param[in] k power of harmonic operation (1: harmonic, 2: biharmonic, etc)
69+
@return Q #V by #V discrete (integrated) k-Laplacian)");
70+
}

tests/test_all.py

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -157,6 +157,13 @@ def test_harmonic():
157157
W = igl.bbw(V,F,b,bc)
158158
W = igl.harmonic(V,F,b,bc,k=1)
159159
W = igl.harmonic(V,F,b,bc,k=2)
160+
161+
def test_harmonic_integrated_from_laplacian_and_mass():
162+
V, F = triangulated_square()
163+
L = igl.cotmatrix(V, F)
164+
M = igl.massmatrix(V, F, igl.MASSMATRIX_TYPE_VORONOI)
165+
Q = igl.harmonic_integrated_from_laplacian_and_mass(L, M, k=1)
166+
Q = igl.harmonic_integrated_from_laplacian_and_mass(L, M, k=2)
160167

161168
def test_tets():
162169
V,F,T = single_tet()

0 commit comments

Comments
 (0)