bug #1330: Cholmod supports double precision only, so let's trigger a static assertion if the scalar type does not match this requirement.

This commit is contained in:
Gael Guennebaud 2016-11-03 10:21:59 +01:00
parent 3e37166d0b
commit 78e93ac1ad
2 changed files with 47 additions and 30 deletions

View File

@ -14,34 +14,40 @@ namespace Eigen {
namespace internal { namespace internal {
template<typename Scalar, typename CholmodType> template<typename Scalar> struct cholmod_configure_matrix;
void cholmod_configure_matrix(CholmodType& mat)
{ template<> struct cholmod_configure_matrix<double> {
if (internal::is_same<Scalar,float>::value) template<typename CholmodType>
{ static void run(CholmodType& mat) {
mat.xtype = CHOLMOD_REAL;
mat.dtype = CHOLMOD_SINGLE;
}
else if (internal::is_same<Scalar,double>::value)
{
mat.xtype = CHOLMOD_REAL; mat.xtype = CHOLMOD_REAL;
mat.dtype = CHOLMOD_DOUBLE; mat.dtype = CHOLMOD_DOUBLE;
} }
else if (internal::is_same<Scalar,std::complex<float> >::value) };
{
mat.xtype = CHOLMOD_COMPLEX; template<> struct cholmod_configure_matrix<std::complex<double> > {
mat.dtype = CHOLMOD_SINGLE; template<typename CholmodType>
} static void run(CholmodType& mat) {
else if (internal::is_same<Scalar,std::complex<double> >::value)
{
mat.xtype = CHOLMOD_COMPLEX; mat.xtype = CHOLMOD_COMPLEX;
mat.dtype = CHOLMOD_DOUBLE; mat.dtype = CHOLMOD_DOUBLE;
} }
else };
{
eigen_assert(false && "Scalar type not supported by CHOLMOD"); // Other scalar types are not yet suppotred by Cholmod
} // template<> struct cholmod_configure_matrix<float> {
} // template<typename CholmodType>
// static void run(CholmodType& mat) {
// mat.xtype = CHOLMOD_REAL;
// mat.dtype = CHOLMOD_SINGLE;
// }
// };
//
// template<> struct cholmod_configure_matrix<std::complex<float> > {
// template<typename CholmodType>
// static void run(CholmodType& mat) {
// mat.xtype = CHOLMOD_COMPLEX;
// mat.dtype = CHOLMOD_SINGLE;
// }
// };
} // namespace internal } // namespace internal
@ -88,7 +94,7 @@ cholmod_sparse viewAsCholmod(SparseMatrix<_Scalar,_Options,_StorageIndex>& mat)
} }
// setup res.xtype // setup res.xtype
internal::cholmod_configure_matrix<_Scalar>(res); internal::cholmod_configure_matrix<_Scalar>::run(res);
res.stype = 0; res.stype = 0;
@ -131,7 +137,7 @@ cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat)
res.x = (void*)(mat.derived().data()); res.x = (void*)(mat.derived().data());
res.z = 0; res.z = 0;
internal::cholmod_configure_matrix<Scalar>(res); internal::cholmod_configure_matrix<Scalar>::run(res);
return res; return res;
} }
@ -180,14 +186,16 @@ class CholmodBase : public SparseSolverBase<Derived>
CholmodBase() CholmodBase()
: m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false)
{ {
m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
cholmod_start(&m_cholmod); cholmod_start(&m_cholmod);
} }
explicit CholmodBase(const MatrixType& matrix) explicit CholmodBase(const MatrixType& matrix)
: m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false)
{ {
m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
cholmod_start(&m_cholmod); cholmod_start(&m_cholmod);
compute(matrix); compute(matrix);
} }
@ -254,7 +262,7 @@ class CholmodBase : public SparseSolverBase<Derived>
eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>()); cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod); cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod);
// If the factorization failed, minor is the column at which it did. On success minor == n. // If the factorization failed, minor is the column at which it did. On success minor == n.
this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue); this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
m_factorizationIsOk = true; m_factorizationIsOk = true;
@ -324,7 +332,7 @@ class CholmodBase : public SparseSolverBase<Derived>
*/ */
Derived& setShift(const RealScalar& offset) Derived& setShift(const RealScalar& offset)
{ {
m_shiftOffset[0] = offset; m_shiftOffset[0] = double(offset);
return derived(); return derived();
} }
@ -386,7 +394,7 @@ class CholmodBase : public SparseSolverBase<Derived>
protected: protected:
mutable cholmod_common m_cholmod; mutable cholmod_common m_cholmod;
cholmod_factor* m_cholmodFactor; cholmod_factor* m_cholmodFactor;
RealScalar m_shiftOffset[2]; double m_shiftOffset[2];
mutable ComputationInfo m_info; mutable ComputationInfo m_info;
int m_factorizationIsOk; int m_factorizationIsOk;
int m_analysisIsOk; int m_analysisIsOk;
@ -410,6 +418,8 @@ class CholmodBase : public SparseSolverBase<Derived>
* *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed. * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
* *
* \warning Only double precision real and complex scalar types are supported by Cholmod.
*
* \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLLT * \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLLT
*/ */
template<typename _MatrixType, int _UpLo = Lower> template<typename _MatrixType, int _UpLo = Lower>
@ -459,6 +469,8 @@ class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimpl
* *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed. * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
* *
* \warning Only double precision real and complex scalar types are supported by Cholmod.
*
* \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLDLT * \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLDLT
*/ */
template<typename _MatrixType, int _UpLo = Lower> template<typename _MatrixType, int _UpLo = Lower>
@ -506,6 +518,8 @@ class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimp
* *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed. * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
* *
* \warning Only double precision real and complex scalar types are supported by Cholmod.
*
* \sa \ref TutorialSparseSolverConcept * \sa \ref TutorialSparseSolverConcept
*/ */
template<typename _MatrixType, int _UpLo = Lower> template<typename _MatrixType, int _UpLo = Lower>
@ -555,6 +569,8 @@ class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSuper
* *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed. * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
* *
* \warning Only double precision real and complex scalar types are supported by Cholmod.
*
* \sa \ref TutorialSparseSolverConcept * \sa \ref TutorialSparseSolverConcept
*/ */
template<typename _MatrixType, int _UpLo = Lower> template<typename _MatrixType, int _UpLo = Lower>

View File

@ -100,7 +100,8 @@
MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY, MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY,
THIS_TYPE_IS_NOT_SUPPORTED, THIS_TYPE_IS_NOT_SUPPORTED,
STORAGE_KIND_MUST_MATCH, STORAGE_KIND_MUST_MATCH,
STORAGE_INDEX_MUST_MATCH STORAGE_INDEX_MUST_MATCH,
CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY
}; };
}; };