Skip to content

Commit be96204

Browse files
authored
Merge pull request #105 from jaganmn/master
Patch integer overflow in 'Rcpp::wrap(<Eigen::Matrix>)'
2 parents 106c26b + a1f7dad commit be96204

File tree

2 files changed

+18
-11
lines changed

2 files changed

+18
-11
lines changed

ChangeLog

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,7 @@
1+
2022-01-15 Mikael Jagan <[email protected]>
2+
3+
* inst/include/RcppEigenWrap.h: Use R_xlen_t for vectors rows + cols
4+
15
2021-12-29 Dirk Eddelbuettel <[email protected]>
26

37
* README.md: Add total downloads badge

inst/include/RcppEigenWrap.h

Lines changed: 14 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -2,7 +2,7 @@
22
//
33
// RcppEigenWrap.h: Rcpp wrap methods for Eigen matrices, vectors and arrays
44
//
5-
// Copyright (C) 2011 - 2012 Douglas Bates, Dirk Eddelbuettel and Romain Francois
5+
// Copyright (C) 2011 - 2022 Douglas Bates, Dirk Eddelbuettel and Romain Francois
66
//
77
// This file is part of RcppEigen.
88
//
@@ -80,16 +80,19 @@ namespace Rcpp{
8080

8181
// for plain dense objects
8282
template <typename T>
83-
SEXP eigen_wrap_plain_dense( const T& obj, Rcpp::traits::true_type ){
84-
typename Eigen::internal::conditional<T::IsRowMajor,
85-
Eigen::Matrix<typename T::Scalar,
86-
T::RowsAtCompileTime,
87-
T::ColsAtCompileTime>,
88-
const T&>::type objCopy(obj);
89-
int m = obj.rows(), n = obj.cols();
90-
R_xlen_t size = static_cast<R_xlen_t>(m) * n;
91-
SEXP ans = PROTECT(::Rcpp::wrap(objCopy.data(), objCopy.data() + size));
92-
if( T::ColsAtCompileTime != 1 ) {
83+
SEXP eigen_wrap_plain_dense( const T& obj, Rcpp::traits::true_type ) {
84+
typename Eigen::internal::conditional<
85+
T::IsRowMajor,
86+
Eigen::Matrix<typename T::Scalar,
87+
T::RowsAtCompileTime,
88+
T::ColsAtCompileTime>,
89+
const T&>::type objCopy(obj);
90+
R_xlen_t m = obj.rows(), n = obj.cols(), size = m * n;
91+
SEXP ans = PROTECT(::Rcpp::wrap(objCopy.data(), objCopy.data() + size));
92+
if ( T::ColsAtCompileTime != 1 ) {
93+
if (m > INT_MAX || n > INT_MAX) {
94+
Rcpp::stop("array dimensions cannot exceed INT_MAX");
95+
}
9396
SEXP dd = PROTECT(::Rf_allocVector(INTSXP, 2));
9497
int *d = INTEGER(dd);
9598
d[0] = m;

0 commit comments

Comments
 (0)