CGAL 5.2.1 - CGAL and Solvers
|
#include <CGAL/Eigen_solver_traits.h>
The class Eigen_solver_traits
provides an interface to the sparse solvers of Eigen.
Eigen version 3.1 (or later) must be available on the system.
EigenSolverT | A sparse solver of Eigen. The default solver is the iterative bi-congugate gradient stabilized solver Eigen::BiCGSTAB for double . |
CGAL::Eigen_sparse_matrix<T>
CGAL::Eigen_sparse_symmetric_matrix<T>
CGAL::Eigen_vector<T>
Instantiation Example
The instantiation of this class assumes an Eigen sparse solver is provided. Here are few examples:
Public Member Functions | |
Eigen_solver_traits () | |
Constructor. | |
bool | linear_solver (const Matrix &A, const Vector &B, Vector &X, NT &D) |
Solve the sparse linear system \( A \times X = B \). More... | |
bool | factor (const Matrix &A, NT &D) |
Factorize the sparse matrix \( A \). More... | |
bool | linear_solver (const Vector &B, Vector &X) |
Solve the sparse linear system \( A \times X = B\), with \( A \) being the matrix provided in factor() . More... | |
bool | linear_solver (const Matrix &B, Vector &X) |
Solve the sparse linear system \( A \times X = B\), with \( A \) being the matrix provided in factor() . More... | |
bool | normal_equation_factor (const Matrix &A) |
Factorize the sparse matrix \( A^t \times A\), where \( A^t \) is the transpose matrix of \( A \). More... | |
bool | normal_equation_solver (const Vector &B, Vector &X) |
Solve the sparse linear system \( A^t \times A \times X = A^t \times B \), with \( A \) being the matrix provided in normal_equation_factor() , and \( A^t \) its transpose matrix. More... | |
bool | normal_equation_solver (const Matrix &A, const Vector &B, Vector &X) |
Equivalent to a call to normal_equation_factor(A) followed by a call to normal_equation_solver(B, X) . More... | |
Protected Attributes | |
const Matrix::EigenType * | m_mat |
boost::shared_ptr< EigenSolverT > | m_solver_sptr |
Types | |
typedef EigenSolverT | Solver |
typedef Scalar | NT |
typedef CGAL::Eigen_vector< NT > | Vector |
typedef Eigen::DenseIndex | Index |
typedef unspecified_type | Matrix |
If T is Eigen::ConjugateGradient<M> or Eigen::SimplicialCholesky<M> , Matrix is CGAL::Eigen_sparse_symmetric_matrix<T> , and CGAL::Eigen_sparse_matrix<T> otherwise. More... | |
Operations | |
EigenSolverT & | solver () |
Returns a reference to the internal Eigen solver. More... | |
typedef unspecified_type CGAL::Eigen_solver_traits< EigenSolverT >::Matrix |
If T
is Eigen::ConjugateGradient<M>
or Eigen::SimplicialCholesky<M>
, Matrix
is CGAL::Eigen_sparse_symmetric_matrix<T>
, and CGAL::Eigen_sparse_matrix<T>
otherwise.
bool CGAL::Eigen_solver_traits< EigenSolverT >::factor | ( | const Matrix & | A, |
NT & | D | ||
) |
Factorize the sparse matrix \( A \).
This factorization is used in linear_solver()
to solve the system for different right-hand side vectors. See linear_solver()
for the description of \( D \).
true
if the factorization is successful and false
otherwise. bool CGAL::Eigen_solver_traits< EigenSolverT >::linear_solver | ( | const Matrix & | A, |
const Vector & | B, | ||
Vector & | X, | ||
NT & | D | ||
) |
Solve the sparse linear system \( A \times X = B \).
Return true
on success. The solution is then \( (1/D) \times X \).
bool CGAL::Eigen_solver_traits< EigenSolverT >::linear_solver | ( | const Vector & | B, |
Vector & | X | ||
) |
Solve the sparse linear system \( A \times X = B\), with \( A \) being the matrix provided in factor()
.
true
if the solver is successful and false
otherwise. bool CGAL::Eigen_solver_traits< EigenSolverT >::linear_solver | ( | const Matrix & | B, |
Vector & | X | ||
) |
Solve the sparse linear system \( A \times X = B\), with \( A \) being the matrix provided in factor()
.
true
if the solver is successful and false
otherwise. bool CGAL::Eigen_solver_traits< EigenSolverT >::normal_equation_factor | ( | const Matrix & | A | ) |
Factorize the sparse matrix \( A^t \times A\), where \( A^t \) is the transpose matrix of \( A \).
This factorization is used in normal_equation_solver()
to solve the system for different right-hand side vectors.
true
if the factorization is successful and false
otherwise. bool CGAL::Eigen_solver_traits< EigenSolverT >::normal_equation_solver | ( | const Vector & | B, |
Vector & | X | ||
) |
Solve the sparse linear system \( A^t \times A \times X = A^t \times B \), with \( A \) being the matrix provided in normal_equation_factor()
, and \( A^t \) its transpose matrix.
true
if the solver is successful and false
otherwise. bool CGAL::Eigen_solver_traits< EigenSolverT >::normal_equation_solver | ( | const Matrix & | A, |
const Vector & | B, | ||
Vector & | X | ||
) |
Equivalent to a call to normal_equation_factor(A)
followed by a call to normal_equation_solver(B, X)
.
EigenSolverT& CGAL::Eigen_solver_traits< EigenSolverT >::solver | ( | ) |
Returns a reference to the internal Eigen solver.
This function can be used for example to set specific parameters of the solver.