From b1eca55328436f9778254c6bbc8852910a0d3d2a Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Tue, 3 Feb 2015 23:46:05 +0100 Subject: Use Ref<> to ensure that both x and b in Ax=b are compatible with Umfpack/SuperLU expectations --- Eigen/src/UmfPackSupport/UmfPackSupport.h | 13 ++++++++++++- 1 file changed, 12 insertions(+), 1 deletion(-) (limited to 'Eigen/src/UmfPackSupport') diff --git a/Eigen/src/UmfPackSupport/UmfPackSupport.h b/Eigen/src/UmfPackSupport/UmfPackSupport.h index b8b216d5e..a2bb75b09 100644 --- a/Eigen/src/UmfPackSupport/UmfPackSupport.h +++ b/Eigen/src/UmfPackSupport/UmfPackSupport.h @@ -403,11 +403,22 @@ bool UmfPackLU::_solve_impl(const MatrixBase &b, MatrixBas eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve"); int errorCode; + Scalar* x_ptr = 0; + Matrix x_tmp; + if(x.innerStride()!=1) + { + x_tmp.resize(x.rows()); + x_ptr = x_tmp.data(); + } for (int j=0; j