Skip to content

Commit c1d915b

Browse files
committed
Split module.cpp
1 parent 63c54cb commit c1d915b

File tree

4 files changed

+43
-33
lines changed

4 files changed

+43
-33
lines changed

CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -82,7 +82,7 @@ target_include_directories(
8282
)
8383
target_link_libraries(nanoeigenpy_headers INTERFACE Eigen3::Eigen)
8484

85-
file(GLOB ${PROJECT_NAME}_SOURCES src/*.cpp)
85+
set(${PROJECT_NAME}_SOURCES src/module.cpp src/solvers.cpp)
8686
nanobind_add_module(nanoeigenpy NB_STATIC NB_SUPPRESS_WARNINGS ${nanoeigenpy_SOURCES} ${nanoeigenpy_HEADERS})
8787
target_link_libraries(nanoeigenpy PRIVATE nanoeigenpy_headers)
8888

src/internal.h

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,9 @@
1+
#include <nanobind/nanobind.h>
2+
#include <Eigen/Core>
3+
4+
namespace nb = nanobind;
5+
6+
using Scalar = double;
7+
static constexpr int Options = Eigen::ColMajor;
8+
using Matrix = Eigen::Matrix<Scalar, -1, -1, Options>;
9+
using Vector = Eigen::Matrix<Scalar, -1, 1>;

src/module.cpp

Lines changed: 2 additions & 32 deletions
Original file line numberDiff line numberDiff line change
@@ -1,21 +1,16 @@
11
/// Copyright 2025 INRIA
22

3-
#include <nanobind/nanobind.h>
43
#include <nanobind/stl/string.h>
54

65
#include "nanoeigenpy/decompositions.hpp"
76
#include "nanoeigenpy/geometry.hpp"
8-
#include "nanoeigenpy/solvers.hpp"
97
#include "nanoeigenpy/utils/is-approx.hpp"
108
#include "nanoeigenpy/constants.hpp"
119

10+
#include "./internal.h"
11+
1212
using namespace nanoeigenpy;
13-
namespace nb = nanobind;
1413

15-
using Scalar = double;
16-
static constexpr int Options = Eigen::ColMajor;
17-
using Matrix = Eigen::Matrix<Scalar, -1, -1, Options>;
18-
using Vector = Eigen::Matrix<Scalar, -1, 1>;
1914
using Quaternion = Eigen::Quaternion<Scalar, Options>;
2015
using SparseMatrix = Eigen::SparseMatrix<Scalar, Options>;
2116

@@ -84,28 +79,3 @@ NB_MODULE(nanoeigenpy, m) {
8479
"Get the set of SIMD instructions used in Eigen when this module was "
8580
"compiled.");
8681
}
87-
88-
void exposeSolvers(nb::module_& m) {
89-
exposeIdentityPreconditioner<Scalar>(m, "IdentityPreconditioner");
90-
exposeDiagonalPreconditioner<Scalar>(m, "DiagonalPreconditioner");
91-
#if EIGEN_VERSION_AT_LEAST(3, 3, 5)
92-
exposeLeastSquareDiagonalPreconditioner<Scalar>(
93-
m, "LeastSquareDiagonalPreconditioner");
94-
#endif
95-
96-
// Solvers
97-
using Eigen::Lower;
98-
using Eigen::Upper;
99-
using ConjugateGradient = Eigen::ConjugateGradient<Matrix, Lower | Upper>;
100-
exposeConjugateGradient<ConjugateGradient>(m, "ConjugateGradient");
101-
using LeastSquaresConjugateGradient = Eigen::LeastSquaresConjugateGradient<
102-
Matrix, Eigen::LeastSquareDiagonalPreconditioner<Scalar>>;
103-
exposeLeastSquaresConjugateGradient<LeastSquaresConjugateGradient>(
104-
m, "LeastSquaresConjugateGradient");
105-
106-
using IdentityConjugateGradient =
107-
Eigen::ConjugateGradient<Matrix, Lower | Upper,
108-
Eigen::IdentityPreconditioner>;
109-
exposeConjugateGradient<IdentityConjugateGradient>(
110-
m, "IdentityConjugateGradient");
111-
}

src/solvers.cpp

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,31 @@
1+
#include "nanoeigenpy/solvers.hpp"
2+
3+
#include "./internal.h"
4+
5+
using namespace nanoeigenpy;
6+
7+
void exposeSolvers(nb::module_& m) {
8+
exposeIdentityPreconditioner<Scalar>(m, "IdentityPreconditioner");
9+
exposeDiagonalPreconditioner<Scalar>(m, "DiagonalPreconditioner");
10+
#if EIGEN_VERSION_AT_LEAST(3, 3, 5)
11+
exposeLeastSquareDiagonalPreconditioner<Scalar>(
12+
m, "LeastSquareDiagonalPreconditioner");
13+
#endif
14+
15+
// Solvers
16+
using Eigen::Lower;
17+
using Eigen::Upper;
18+
using ConjugateGradient = Eigen::ConjugateGradient<Matrix, Lower | Upper>;
19+
exposeConjugateGradient<ConjugateGradient>(m, "ConjugateGradient");
20+
21+
using LeastSquaresConjugateGradient = Eigen::LeastSquaresConjugateGradient<
22+
Matrix, Eigen::LeastSquareDiagonalPreconditioner<Scalar>>;
23+
exposeLeastSquaresConjugateGradient<LeastSquaresConjugateGradient>(
24+
m, "LeastSquaresConjugateGradient");
25+
26+
using IdentityConjugateGradient =
27+
Eigen::ConjugateGradient<Matrix, Lower | Upper,
28+
Eigen::IdentityPreconditioner>;
29+
exposeConjugateGradient<IdentityConjugateGradient>(
30+
m, "IdentityConjugateGradient");
31+
}

0 commit comments

Comments
 (0)