aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseCholesky
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-10-07 18:29:28 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-10-07 18:29:28 +0200
commit57413492948d4c9dd6572f5bde6720b80122054a (patch)
tree22c8b070d5f929c32428e562505a27f704d978eb /Eigen/src/SparseCholesky
parent118b1113d9a4fa1a263597cdd699e237e5f2b4ac (diff)
bug #882: fix various const-correctness issues with *View classes.
Diffstat (limited to 'Eigen/src/SparseCholesky')
-rw-r--r--Eigen/src/SparseCholesky/SimplicialCholesky.h8
1 files changed, 4 insertions, 4 deletions
diff --git a/Eigen/src/SparseCholesky/SimplicialCholesky.h b/Eigen/src/SparseCholesky/SimplicialCholesky.h
index 0e8fa6628..918a34e13 100644
--- a/Eigen/src/SparseCholesky/SimplicialCholesky.h
+++ b/Eigen/src/SparseCholesky/SimplicialCholesky.h
@@ -237,8 +237,8 @@ template<typename _MatrixType, int _UpLo, typename _Ordering> struct traits<Simp
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::Index Index;
typedef SparseMatrix<Scalar, ColMajor, Index> CholMatrixType;
- typedef TriangularView<CholMatrixType, Eigen::Lower> MatrixL;
- typedef TriangularView<typename CholMatrixType::AdjointReturnType, Eigen::Upper> MatrixU;
+ typedef TriangularView<const CholMatrixType, Eigen::Lower> MatrixL;
+ typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::Upper> MatrixU;
static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
};
@@ -251,8 +251,8 @@ template<typename _MatrixType,int _UpLo, typename _Ordering> struct traits<Simpl
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::Index Index;
typedef SparseMatrix<Scalar, ColMajor, Index> CholMatrixType;
- typedef TriangularView<CholMatrixType, Eigen::UnitLower> MatrixL;
- typedef TriangularView<typename CholMatrixType::AdjointReturnType, Eigen::UnitUpper> MatrixU;
+ typedef TriangularView<const CholMatrixType, Eigen::UnitLower> MatrixL;
+ typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::UnitUpper> MatrixU;
static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
};