CGAL 5.5.2 - CGAL and Solvers
CGAL::Eigen_solver_traits< EigenSolverT > Class Template Reference

#include <CGAL/Eigen_solver_traits.h>

Definition

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
class CGAL::Eigen_solver_traits< EigenSolverT >

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.

Is Model Of:
SparseLinearAlgebraWithFactorTraits_d and NormalEquationSparseLinearAlgebraTraits_d
Template Parameters
EigenSolverTA sparse solver of Eigen. The default solver is the iterative bi-conjugate gradient stabilized solver Eigen::BiCGSTAB for double.
See also
CGAL::Eigen_sparse_matrix<T>
CGAL::Eigen_sparse_symmetric_matrix<T>
CGAL::Eigen_vector<T>
http://eigen.tuxfamily.org/index.php?title=Main_Page

Instantiation Example

The instantiation of this class assumes an Eigen sparse solver is provided. Here are few examples:

//iterative general solver
//iterative symmetric solver
//direct symmetric solver
Examples:
Solver_interface/sparse_solvers.cpp.

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
 
std::shared_ptr< EigenSolverT > m_solver_sptr
 

Types

typedef EigenSolverT Solver
 
typedef Scalar NT
 
typedef CGAL::Eigen_vector< NTVector
 
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...
 

Member Typedef Documentation

◆ Matrix

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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.

Member Function Documentation

◆ factor()

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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 \).

Returns
true if the factorization is successful and false otherwise.

◆ linear_solver() [1/3]

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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 \).

Precondition
A.row_dimension() == B.dimension().
A.column_dimension() == X.dimension().

◆ linear_solver() [2/3]

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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().

Returns
true if the solver is successful and false otherwise.

◆ linear_solver() [3/3]

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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().

Returns
true if the solver is successful and false otherwise.

◆ normal_equation_factor()

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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.

Returns
true if the factorization is successful and false otherwise.

◆ normal_equation_solver() [1/2]

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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.

Returns
true if the solver is successful and false otherwise.

◆ normal_equation_solver() [2/2]

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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) .

◆ solver()

template<class EigenSolverT = Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType>>
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.