Include Lapack svd in Solver interface (another model for SvdTraits)

This commit is contained in:
Simon Giraudot 2015-08-26 14:48:11 +02:00
parent 7e50599ffa
commit 8a46d3cd5e
6 changed files with 28 additions and 11 deletions

View File

@ -508,6 +508,15 @@ and the \ref PkgRidges_3 packages.
The \sc{Eigen} web site is <A HREF="http://eigen.tuxfamily.org">`http://eigen.tuxfamily.org`</A>.
\subsection thirdpartyLapack Lapack
\sc{Lapack} is a `Fortran` library for linear algebra.
In \cgal, \sc{Lapack} is an alternative to \sc{Eigen} for singular value decomposition for the \ref
kgJet_fitting_3 and the \ref PkgRidges_3 packages.
The \sc{Lapack} web site is <A HREF="http://www.netlib.org/lapack/">`http://www.netlib.org/lapack/`</A>.
\subsection thirdpartylibQGLViewer libQGLViewer
libQGLViewer is a 3D widget based on \sc{Qt} 4's `QGLWidget`. In case of \sc{Qt}5 used, libQGLViewer needs to be recompiled with the proper \sc{Qt}5 version.

View File

@ -42,15 +42,15 @@ class:
The concept `SvdTraits` defines an interface for solving in the least
square sense a linear system with a singular value decomposition. The
field type is `double`.
field type is `double`. We provide two models for this concept:
The class `Eigen_svd` provides an implementation of `SvdTraits` using
the \ref thirdpartyEigen library.
- `Eigen_svd` uses the \ref thirdpartyEigen library.
- `Lapack_svd` uses the \ref thirdpartyLapack library.
Here is a simple example that shows how to handle matrices, vectors
and this solver:
\cgalExample{Solver_interface/eigen_singular_value_decomposition.cpp}
\cgalExample{Solver_interface/singular_value_decomposition.cpp}

View File

@ -1,5 +1,5 @@
/*!
\example Solver_interface/diagonalize_matrix.cpp
\example Solver_interface/eigen_singular_value_decomposition.cpp
\example Solver_interface/singular_value_decomposition.cpp
\example Solver_interface/eigen_sparse_solvers.cpp
*/

View File

@ -24,7 +24,7 @@ if ( CGAL_FOUND )
include_directories (BEFORE "../include")
create_single_source_cgal_program( "eigen_singular_value_decomposition.cpp" )
create_single_source_cgal_program( "singular_value_decomposition.cpp" )
create_single_source_cgal_program( "diagonalize_matrix.cpp" )
create_single_source_cgal_program( "eigen_sparse_solvers.cpp" )

View File

@ -1,11 +1,20 @@
#include <CGAL/Simple_cartesian.h>
#ifdef CGAL_EIGEN3_ENABLED
#include <CGAL/Eigen_matrix.h>
#include <CGAL/Eigen_vector.h>
#include <CGAL/Eigen_svd.h>
typedef CGAL::Eigen_svd Svd;
#else
#ifdef CGAL_LAPACK_ENABLED
#include <CGAL/Lapack_svd.h>
typedef CGAL::Lapack_svd Svd;
#endif
#endif
typedef CGAL::Eigen_svd::FT FT;
typedef CGAL::Eigen_svd::Vector Eigen_vector;
typedef CGAL::Eigen_svd::Matrix Eigen_matrix;
typedef Svd::FT FT;
typedef Svd::Vector Eigen_vector;
typedef Svd::Matrix Eigen_matrix;
int main(void)
{
@ -14,7 +23,6 @@ int main(void)
Eigen_vector B (degree);
Eigen_matrix M (degree, degree);
// Fill B and M with random numbers
for (std::size_t i = 0; i < degree; ++ i)
{
@ -24,7 +32,7 @@ int main(void)
}
// Solve MX=B
std::cout << CGAL::Eigen_svd::solve(M, B) << std::endl;
std::cout << Svd::solve(M, B) << std::endl;
// Print result
std::cout << "Solution of SVD = [ ";