Skip to content
Snippets Groups Projects
Verified Commit df49d6ee authored by Justin Carpentier's avatar Justin Carpentier
Browse files

test/sparse: add minimal test

parent 77c8684c
No related branches found
No related tags found
No related merge requests found
Pipeline #35628 passed with warnings
# #
# Copyright (c) 2014-2019 CNRS Copyright (c) 2018-2023 INRIA # Copyright (c) 2014-2019 CNRS Copyright (c) 2018-2024 INRIA
# #
macro(ADD_LIB_UNIT_TEST test) macro(ADD_LIB_UNIT_TEST test)
...@@ -25,6 +25,7 @@ macro(ADD_LIB_UNIT_TEST test) ...@@ -25,6 +25,7 @@ macro(ADD_LIB_UNIT_TEST test)
endmacro(ADD_LIB_UNIT_TEST) endmacro(ADD_LIB_UNIT_TEST)
add_lib_unit_test(matrix) add_lib_unit_test(matrix)
add_lib_unit_test(sparse_matrix)
add_lib_unit_test(tensor) add_lib_unit_test(tensor)
add_lib_unit_test(geometry) add_lib_unit_test(geometry)
add_lib_unit_test(complex) add_lib_unit_test(complex)
......
from __future__ import print_function
import numpy as np
import sparse_matrix
m = sparse_matrix.emptyMatrix()
assert m.shape == (0, 0)
v = sparse_matrix.emptyVector()
assert v.shape == (0, 0)
m = sparse_matrix.matrix1x1(2)
assert m.toarray() == np.array([2])
v = sparse_matrix.vector1x1(2)
assert v.toarray() == np.array([2])
size = 100
diag_values = np.random.rand(100)
diag_mat = sparse_matrix.diagonal(diag_values)
assert (diag_mat.toarray() == np.diag(diag_values)).all()
diag_mat_copy = sparse_matrix.copy(diag_mat)
assert (diag_mat_copy != diag_mat).nnz == 0
/*
* Copyright 2024 CNRS INRIA
*/
#include <iostream>
#include "eigenpy/eigenpy.hpp"
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> vector1x1(const Scalar& value) {
typedef Eigen::SparseMatrix<Scalar, Options> ReturnType;
ReturnType mat(1, 1);
mat.coeffRef(0, 0) = value;
mat.makeCompressed();
return mat;
}
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> matrix1x1(const Scalar& value) {
typedef Eigen::SparseMatrix<Scalar, Options> ReturnType;
ReturnType mat(1, 1);
mat.coeffRef(0, 0) = value;
mat.makeCompressed();
return mat;
}
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> diagonal(
const Eigen::Ref<Eigen::Matrix<Scalar, Eigen::Dynamic, 1> >& diag_values) {
typedef Eigen::SparseMatrix<Scalar, Options> ReturnType;
ReturnType mat(diag_values.size(), diag_values.size());
for (Eigen::Index k = 0; k < diag_values.size(); ++k)
mat.coeffRef(k, k) = diag_values[k];
mat.makeCompressed();
return mat;
}
template <typename Scalar, int Options>
void matrix1x1_input(const Eigen::Matrix<Scalar, 1, 1>& mat) {
std::cout << mat << std::endl;
}
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> emptyVector() {
return Eigen::SparseMatrix<Scalar, Options>();
}
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> emptyMatrix() {
return Eigen::SparseMatrix<Scalar, Options>();
}
template <typename Scalar, int Options>
void print(const Eigen::SparseMatrix<Scalar, Options>& mat) {
std::cout << mat << std::endl;
}
template <typename Scalar, int Options>
Eigen::SparseMatrix<Scalar, Options> copy(
const Eigen::SparseMatrix<Scalar, Options>& mat) {
return mat;
}
BOOST_PYTHON_MODULE(sparse_matrix) {
using namespace Eigen;
namespace bp = boost::python;
eigenpy::enableEigenPy();
typedef Eigen::SparseMatrix<double> SparseMatrixD;
eigenpy::EigenToPyConverter<SparseMatrixD>::registration();
eigenpy::EigenFromPyConverter<SparseMatrixD>::registration();
bp::def("vector1x1", vector1x1<double, Eigen::ColMajor>);
bp::def("matrix1x1", matrix1x1<double, Eigen::ColMajor>);
bp::def("matrix1x1", matrix1x1_input<double, Eigen::ColMajor>);
bp::def("matrix1x1_int", matrix1x1_input<int, Eigen::ColMajor>);
bp::def("print", print<double, Eigen::ColMajor>);
bp::def("copy", copy<double, Eigen::ColMajor>);
bp::def("diagonal", diagonal<double, Eigen::ColMajor>);
bp::def("emptyVector", emptyVector<double, Eigen::ColMajor>);
bp::def("emptyMatrix", emptyMatrix<double, Eigen::ColMajor>);
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment