Skip to content
Snippets Groups Projects
Commit 824a52b3 authored by jcarpent's avatar jcarpent
Browse files

[Solvers] Expose iterative solvers of Eigen

parent c81dd85c
No related branches found
No related tags found
No related merge requests found
#
# Copyright (c) 2015-2016 LAAS-CNRS
# Copyright (c) 2015-2017 LAAS-CNRS
#
# This file is part of eigenpy.
# eigenpy is free software: you can redistribute it and/or
......@@ -90,8 +90,13 @@ SET(HEADERS
registration.hpp
angle-axis.hpp
quaternion.hpp
solvers/solvers.hpp
solvers/IterativeSolverBase.hpp
solvers/ConjugateGradient.hpp
solvers/SparseSolverBase.hpp
)
MAKE_DIRECTORY("${${PROJECT_NAME}_BINARY_DIR}/include/eigenpy")
MAKE_DIRECTORY("${${PROJECT_NAME}_BINARY_DIR}/include/eigenpy/solvers")
INCLUDE_DIRECTORIES(${${PROJECT_NAME}_BINARY_DIR}/include/eigenpy)
#FOREACH(header ${HEADERS})
......@@ -131,6 +136,7 @@ SET(${PROJECT_NAME}_SOURCES
src/eigenpy.cpp
src/angle-axis.cpp
src/quaternion.cpp
src/solvers/solvers.cpp
)
ADD_LIBRARY(${PROJECT_NAME} SHARED ${${PROJECT_NAME}_SOURCES} ${${PROJECT_NAME}_HEADERS})
......
......@@ -17,6 +17,7 @@
#include "eigenpy/eigenpy.hpp"
#include "eigenpy/geometry.hpp"
#include "eigenpy/solvers/solvers.hpp"
#include <iostream>
......@@ -26,5 +27,6 @@ BOOST_PYTHON_MODULE(eigenpy)
eigenpy::enableEigenPy();
eigenpy::exposeAngleAxis();
eigenpy::exposeQuaternion();
eigenpy::exposeSolvers();
}
/*
* Copyright 2017, Justin Carpentier, LAAS-CNRS
*
* This file is part of eigenpy.
* eigenpy is free software: you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation, either version 3 of
* the License, or (at your option) any later version.
* eigenpy is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied warranty
* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details. You should
* have received a copy of the GNU Lesser General Public License along
* with eigenpy. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef __eigenpy_conjugate_gradient_hpp__
#define __eigenpy_conjugate_gradient_hpp__
#include <boost/python.hpp>
#include <Eigen/IterativeLinearSolvers>
#include "eigenpy/solvers/IterativeSolverBase.hpp"
namespace eigenpy
{
namespace bp = boost::python;
template<typename ConjugateGradient>
struct ConjugateGradientVisitor
: public boost::python::def_visitor< ConjugateGradientVisitor<ConjugateGradient> >
{
typedef typename ConjugateGradient::MatrixType MatrixType;
template<class PyClass>
void visit(PyClass& cl) const
{
cl
.def(bp::init<>("Default constructor"))
.def(bp::init<MatrixType>(bp::arg("A"),"Initialize the solver with matrix A for further Ax=b solving.\n"
"This constructor is a shortcut for the default constructor followed by a call to compute()."))
;
}
static void expose()
{
bp::class_<ConjugateGradient,boost::noncopyable>("ConjugateGradient",
bp::no_init)
.def(IterativeSolverVisitor<ConjugateGradient>())
.def(ConjugateGradientVisitor<ConjugateGradient>())
;
}
};
} // namespace eigenpy
#endif // ifndef __eigenpy_conjugate_gradient_hpp__
/*
* Copyright 2017, Justin Carpentier, LAAS-CNRS
*
* This file is part of eigenpy.
* eigenpy is free software: you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation, either version 3 of
* the License, or (at your option) any later version.
* eigenpy is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied warranty
* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details. You should
* have received a copy of the GNU Lesser General Public License along
* with eigenpy. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef __eigenpy_iterative_solver_base_hpp__
#define __eigenpy_iterative_solver_base_hpp__
#include <boost/python.hpp>
#include <Eigen/Core>
#include "eigenpy/solvers/SparseSolverBase.hpp"
namespace eigenpy
{
namespace bp = boost::python;
template<typename IterativeSolver>
struct IterativeSolverVisitor
: public boost::python::def_visitor< IterativeSolverVisitor<IterativeSolver> >
{
typedef typename IterativeSolver::MatrixType MatrixType;
typedef Eigen::VectorXd VectorType;
template<class PyClass>
void visit(PyClass& cl) const
{
typedef IterativeSolver IS;
cl
.def(SparseSolverVisitor<IS>())
.def("error",&IS::error,"Returns the tolerance error reached during the last solve.\n"
"It is a close approximation of the true relative residual error |Ax-b|/|b|.")
.def("info",&IS::info,"Returns success if the iterations converged, and NoConvergence otherwise.")
.def("iterations",&IS::iterations,"Returns the number of iterations performed during the last solve.")
.def("maxIterations",&IS::maxIterations,"Returns the max number of iterations.\n"
"It is either the value setted by setMaxIterations or, by default, twice the number of columns of the matrix.")
.def("setMaxIterations",&IS::setMaxIterations,"Sets the max number of iterations.\n"
"Default is twice the number of columns of the matrix.",
bp::return_value_policy<bp::reference_existing_object>())
.def("tolerance",&IS::tolerance,"Returns he tolerance threshold used by the stopping criteria.")
.def("setTolerance",&IS::setTolerance,"Sets the tolerance threshold used by the stopping criteria.\n"
"This value is used as an upper bound to the relative residual error: |Ax-b|/|b|. The default value is the machine precision.",
bp::return_value_policy<bp::reference_existing_object>())
.def("analyzePattern",&analyzePattern,bp::arg("A"),"Initializes the iterative solver for the sparsity pattern of the matrix A for further solving Ax=b problems.\n"
"Currently, this function mostly calls analyzePattern on the preconditioner.\n"
"In the future we might, for instance, implement column reordering for faster matrix vector products.",
bp::return_value_policy<bp::reference_existing_object>())
.def("factorize",&factorize,bp::arg("A"),"Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.\n"
"Currently, this function mostly calls factorize on the preconditioner.",
bp::return_value_policy<bp::reference_existing_object>())
.def("compute",&compute,bp::arg("A"),"Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.\n"
"Currently, this function mostly calls factorize on the preconditioner.\n"
"In the future we might, for instance, implement column reordering for faster matrix vector products.",
bp::return_value_policy<bp::reference_existing_object>())
.def("solveWithGuess",&solveWithGuess,bp::args("b","x0"),
"Returns the solution x of Ax = b using the current decomposition of A and x0 as an initial solution.")
;
}
private:
static IterativeSolver & factorize(IterativeSolver & self, const MatrixType & m)
{
return self.factorize(m);
}
static IterativeSolver & compute(IterativeSolver & self, const MatrixType & m)
{
return self.compute(m);
}
static IterativeSolver & analyzePattern(IterativeSolver & self, const MatrixType & m)
{
return self.analyzePattern(m);
}
static VectorType solveWithGuess(IterativeSolver & self, const Eigen::VectorXd & b, const Eigen::VectorXd & x0)
{
return self.solveWithGuess(b,x0);
}
};
} // namespace eigenpy
#endif // ifndef __eigenpy_iterative_solver_base_hpp__
/*
* Copyright 2017, Justin Carpentier, LAAS-CNRS
*
* This file is part of eigenpy.
* eigenpy is free software: you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation, either version 3 of
* the License, or (at your option) any later version.
* eigenpy is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied warranty
* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details. You should
* have received a copy of the GNU Lesser General Public License along
* with eigenpy. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef __eigenpy_sparse_solver_base_hpp__
#define __eigenpy_sparse_solver_base_hpp__
#include <boost/python.hpp>
#include <Eigen/Core>
namespace eigenpy
{
namespace bp = boost::python;
template<typename SparseSolver>
struct SparseSolverVisitor
: public bp::def_visitor< SparseSolverVisitor<SparseSolver> >
{
typedef Eigen::VectorXd VectorType;
template<class PyClass>
void visit(PyClass& cl) const
{
cl
.def("solve",&solve,bp::arg("b"),
"Returns the solution x of Ax = b using the current decomposition of A.")
;
}
private:
static VectorType solve(SparseSolver & self, const VectorType & b)
{
return self.solve(b);
}
};
} // namespace eigenpy
#endif // ifndef __eigenpy_sparse_solver_base_hpp__
/*
* Copyright 2017, Justin Carpentier, LAAS-CNRS
*
* This file is part of eigenpy.
* eigenpy is free software: you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation, either version 3 of
* the License, or (at your option) any later version.
* eigenpy is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied warranty
* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details. You should
* have received a copy of the GNU Lesser General Public License along
* with eigenpy. If not, see <http://www.gnu.org/licenses/>.
*/
#include "eigenpy/solvers/solvers.hpp"
#include "eigenpy/solvers/ConjugateGradient.hpp"
namespace eigenpy
{
void exposeSolvers()
{
using namespace Eigen;
ConjugateGradientVisitor< ConjugateGradient<MatrixXd,Lower|Upper> >::expose();
boost::python::enum_<Eigen::ComputationInfo>("ComputationInfo")
.value("Success",Eigen::Success)
.value("NumericalIssue",Eigen::NumericalIssue)
.value("NoConvergence",Eigen::NoConvergence)
.value("InvalidInput",Eigen::InvalidInput)
;
}
} // namespace eigenpy
/*
* Copyright 2017, Justin Carpentier, LAAS-CNRS
*
* This file is part of eigenpy.
* eigenpy is free software: you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation, either version 3 of
* the License, or (at your option) any later version.
* eigenpy is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied warranty
* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details. You should
* have received a copy of the GNU Lesser General Public License along
* with eigenpy. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef __eigenpy_solvers_hpp__
#define __eigenpy_solvers_hpp__
namespace eigenpy
{
void exposeSolvers();
} // namespace eigenpy
#endif // define __eigenpy_solvers_hpp__
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