Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Header clean-up to fix lme4/lme4#745 #131

Merged
merged 4 commits into from
Nov 1, 2023
Merged
Show file tree
Hide file tree
Changes from 2 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 2 additions & 2 deletions DESCRIPTION
Original file line number Diff line number Diff line change
Expand Up @@ -24,7 +24,7 @@ License: GPL (>= 2) | file LICENSE
LazyLoad: yes
Depends: R (>= 3.6.0)
LinkingTo: Rcpp
Imports: Matrix (>= 1.1-0), Rcpp (>= 0.11.0), stats, utils
Suggests: inline, tinytest, pkgKitten, microbenchmark
Imports: Rcpp (>= 0.11.0), stats, utils
Suggests: Matrix, inline, microbenchmark, pkgKitten, tinytest
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Just a side comment and not a request to change: When I modify code (particular "other people's code" I try to be minimal. Adding Matrix was needed on the Suggests:, doing it at the front is fine, reordering the four remaining words just puts cognitive load here that we do not need. (Again, no need to change, just a comment.)

URL: https://github.com/RcppCore/RcppEigen, https://dirk.eddelbuettel.com/code/rcpp.eigen.html
BugReports: https://github.com/RcppCore/RcppEigen/issues
1 change: 0 additions & 1 deletion NAMESPACE
Original file line number Diff line number Diff line change
@@ -1,6 +1,5 @@
useDynLib("RcppEigen", .registration=TRUE)

importClassesFrom("Matrix", "dgCMatrix", "dgeMatrix", "dsCMatrix", "dtCMatrix")
importFrom("Rcpp", "evalCpp")
importFrom("utils", "packageDescription", "package.skeleton")
importFrom("stats", "model.frame", "model.matrix", "model.response", "fitted", "coef", "printCoefmat", "pt")
Expand Down
4 changes: 0 additions & 4 deletions inst/include/Eigen/CholmodSupport
Original file line number Diff line number Diff line change
Expand Up @@ -12,10 +12,6 @@

#include "src/Core/util/DisableStupidWarnings.h"

extern "C" {
#include <RcppEigenCholmod.h>
}

/** \ingroup Support_modules
* \defgroup CholmodSupport_Module CholmodSupport module
*
Expand Down
22 changes: 11 additions & 11 deletions inst/include/Eigen/src/CholmodSupport/CholmodSupport.h
Original file line number Diff line number Diff line change
Expand Up @@ -195,23 +195,23 @@ class CholmodBase : public SparseSolverBase<Derived>
{
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);
R_MATRIX_CHOLMOD(start)(&m_cholmod);
}

explicit CholmodBase(const MatrixType& matrix)
: m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false)
{
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);
R_MATRIX_CHOLMOD(start)(&m_cholmod);
compute(matrix);
}

~CholmodBase()
{
if(m_cholmodFactor)
cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
cholmod_finish(&m_cholmod);
R_MATRIX_CHOLMOD(free_factor)(&m_cholmodFactor, &m_cholmod);
R_MATRIX_CHOLMOD(finish)(&m_cholmod);
}

inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
Expand Down Expand Up @@ -246,11 +246,11 @@ class CholmodBase : public SparseSolverBase<Derived>
{
if(m_cholmodFactor)
{
cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
R_MATRIX_CHOLMOD(free_factor)(&m_cholmodFactor, &m_cholmod);
m_cholmodFactor = 0;
}
cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
m_cholmodFactor = R_MATRIX_CHOLMOD(analyze)(&A, &m_cholmod);

this->m_isInitialized = true;
this->m_info = Success;
Expand All @@ -268,7 +268,7 @@ class CholmodBase : public SparseSolverBase<Derived>
{
eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod);
R_MATRIX_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.
this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
Expand All @@ -293,15 +293,15 @@ class CholmodBase : public SparseSolverBase<Derived>
Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b.derived());

cholmod_dense b_cd = viewAsCholmod(b_ref);
cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
cholmod_dense* x_cd = R_MATRIX_CHOLMOD(solve)(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
if(!x_cd)
{
this->m_info = NumericalIssue;
return;
}
// TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
cholmod_free_dense(&x_cd, &m_cholmod);
R_MATRIX_CHOLMOD(free_dense)(&x_cd, &m_cholmod);
}

/** \internal */
Expand All @@ -316,15 +316,15 @@ class CholmodBase : public SparseSolverBase<Derived>
// note: cs stands for Cholmod Sparse
Ref<SparseMatrix<typename RhsDerived::Scalar,ColMajor,typename RhsDerived::StorageIndex> > b_ref(b.const_cast_derived());
cholmod_sparse b_cs = viewAsCholmod(b_ref);
cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
cholmod_sparse* x_cs = R_MATRIX_CHOLMOD(spsolve)(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
if(!x_cs)
{
this->m_info = NumericalIssue;
return;
}
// TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
dest.derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
cholmod_free_sparse(&x_cs, &m_cholmod);
R_MATRIX_CHOLMOD(free_sparse)(&x_cs, &m_cholmod);
}
#endif // EIGEN_PARSED_BY_DOXYGEN

Expand Down
Loading