34 #include <visp3/core/vpConfig.h>
36 #include <visp3/core/vpColVector.h>
37 #include <visp3/core/vpMath.h>
38 #include <visp3/core/vpMatrix.h>
41 #include <visp3/core/vpException.h>
42 #include <visp3/core/vpMatrixException.h>
44 #if defined(VISP_HAVE_OPENCV)
45 #include <opencv2/core/core.hpp>
48 #ifdef VISP_HAVE_LAPACK
50 #include <gsl/gsl_linalg.h>
54 typedef MKL_INT integer;
55 #elif !defined(VISP_HAVE_GSL)
56 #ifdef VISP_HAVE_LAPACK_BUILT_IN
57 typedef long int integer;
61 extern "C" void dpotrf_(
char *uplo, integer *n,
double *a, integer *lda, integer *info);
62 extern "C" int dpotri_(
char *uplo, integer *n,
double *a, integer *lda, integer *info);
66 #if defined(VISP_HAVE_EIGEN3)
67 #include <Eigen/Dense>
114 #if defined(VISP_HAVE_LAPACK)
116 #elif defined(VISP_HAVE_OPENCV)
123 #if defined(VISP_HAVE_LAPACK)
166 #if defined(VISP_HAVE_GSL)
178 #if (GSL_MAJOR_VERSION >= 2 && GSL_MINOR_VERSION >= 3)
179 gsl_linalg_cholesky_decomp1(&
cholesky);
181 gsl_linalg_cholesky_decomp(&
cholesky);
183 gsl_linalg_cholesky_invert(&
cholesky);
193 integer rowNum_ = (integer)this->
getRows();
194 integer lda = (integer)rowNum_;
198 dpotrf_((
char *)
"L", &rowNum_, A.
data, &lda, &info);
201 std::stringstream errMsg;
202 errMsg <<
"Cannot inverse by Cholesky with Lapack: error "<< info <<
" in dpotrf_()";
206 dpotri_((
char *)
"L", &rowNum_, A.
data, &lda, &info);
208 std::cout <<
"cholesky:dpotri_:error" << std::endl;
212 for (
unsigned int i = 0; i < A.
getRows(); ++i)
213 for (
unsigned int j = 0; j < A.
getCols(); ++j)
231 #if defined(VISP_HAVE_GSL)
240 #if (GSL_MAJOR_VERSION >= 2 && GSL_MINOR_VERSION >= 3)
241 gsl_linalg_cholesky_decomp1(&
cholesky);
243 gsl_linalg_cholesky_decomp(&
cholesky);
251 integer rowNum_ = (integer)this->
getRows();
252 integer lda = (integer)rowNum_;
254 dpotrf_((
char *)
"L", &rowNum_, L.data, &lda, &info);
257 std::stringstream errMsg;
258 errMsg <<
"The " << -info <<
"th argument has an illegal value" << std::endl;
262 std::stringstream errMsg;
263 errMsg <<
"The leading minor of order" << info <<
"is not positive definite." << std::endl;
268 unsigned int nbRows = this->
getRows();
269 unsigned int nbCols = this->
getCols();
270 for (
unsigned int r = 0; r < nbRows; ++r) {
271 for (
unsigned int c = r + 1; c < nbCols; ++c) {
279 #if defined(VISP_HAVE_OPENCV)
328 cv::Mat Minv = M.inv(cv::DECOMP_CHOLESKY);
331 memcpy(A.
data, Minv.data, (
size_t)(8 * Minv.rows * Minv.cols));
345 memcpy(M.data, this->data, (
size_t)(8 * M.rows * M.cols));
346 std::size_t bytes_per_row =
sizeof(double)*
rowNum;
347 bool result = cv::Cholesky(M.ptr<
double>(), bytes_per_row,
rowNum,
nullptr, bytes_per_row,
rowNum);
352 memcpy(L.data, M.data, (
size_t)(8 * M.rows * M.cols));
354 unsigned int nbRows = this->
getRows();
355 unsigned int nbCols = this->
getCols();
356 for (
unsigned int r = 0; r < nbRows; ++r) {
357 for (
unsigned int c = r + 1; c < nbCols; ++c) {
365 #if defined(VISP_HAVE_EIGEN3)
374 unsigned int nbRows = this->
getRows();
375 unsigned int nbCols = this->
getCols();
376 Eigen::MatrixXd A(nbRows, nbCols);
377 for (
unsigned int r = 0; r < nbRows; ++r) {
378 for (
unsigned int c = 0; c < nbCols; ++c) {
379 A(r, c) = (*this)[r][c];
382 Eigen::MatrixXd L = A.llt().matrixL();
383 vpMatrix Lvisp(
static_cast<unsigned int>(L.rows()),
static_cast<unsigned int>(L.cols()), 0.);
384 for (
unsigned int r = 0; r < nbRows; ++r) {
385 for (
unsigned int c = 0; c <= r; ++c) {
386 Lvisp[r][c] = L(r, c);
unsigned int getCols() const
Type * data
Address of the first element of the data array.
unsigned int rowNum
Number of rows in the array.
unsigned int getRows() const
unsigned int colNum
Number of columns in the array.
error that can be emitted by ViSP classes.
@ badValue
Used to indicate that a value is not in the allowed range.
error that can be emitted by the vpMatrix class and its derivatives
@ forbiddenOperatorError
Forbidden operation.
@ matrixError
Matrix operation error.
Implementation of a matrix and operations on matrices.
vpMatrix inverseByCholesky() const
vpMatrix cholesky() const
vpMatrix inverseByCholeskyLapack() const
vpMatrix choleskyByOpenCV() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using OpenCV library.
vpMatrix choleskyByEigen3() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Eigen3 library.
vpMatrix inverseByCholeskyOpenCV() const
vpMatrix choleskyByLapack() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Lapack library.