[eigen] Inlining (or lack thereof) on msvc 2013

[ Thread Index | Date Index | More lists.tuxfamily.org/eigen Archives ]


Hello,

we are observing some inlining issues with MSVC from Visual Studio 2013 (Update 5). 
On many occasions the compiler does not inline any of the Eigen function calls.

This seems to be somehow related to the size of the function that is using Eigen.

I have attached an example (eigen_inline.cpp) where the compiler does not inline
any of the Eigen operations (see eigen_inline.asm, search for ENDP ; compute):
....
	call	?functor@?$CwiseUnaryOp@U?$scalar_multiple_op@N@internal@Eigen@@$$CBV?$Block@$$CBV?$Matrix@N$02$00$0A@$02$00@Eigen@@$01$00$0A@@3@@Eigen@@QEBAAEBU?$scalar_multiple_op@N@internal@2@XZ ; Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<double>,Eigen::Block<Eigen::Matrix<double,3,1,0,3,1> const ,2,1,0> const >::functor
.....
That of course results in a considerable performance hit.

However if we remove one single line of source code from the sample function, 
or alternatively replace "inline" with "EIGEN_STRONG_INLINE" in all of Eigen, 
then all calls are properly inlined (see eigen_strong_inline.asm):
....
	movsdx	QWORD PTR tmp2$[rsp], xmm6
	movsdx	xmm1, QWORD PTR [r9+24]
	movsdx	xmm0, QWORD PTR [r9+32]
	mulsd	xmm1, xmm4
	mulsd	xmm0, xmm3
	subsd	xmm0, xmm1

All code was compiled with:
  cl.exe -c /O2 /FA /DNDEBUG -I..\Eigen eigen_inline.cpp

Do you have any insights into this kind of behavior? A less invasive solution that does not include
patching all of Eigen would be preferable...

Best
Simon
#include <Eigen/Core>

using Matrix2x3d = Eigen::Matrix<double, 2, 3>;

Matrix2x3d compute(const Eigen::VectorXd& m, 
                   const Eigen::Vector3d& M, 
				   const Eigen::Matrix3d& R, 
				   Eigen::Vector3d& P)
{
		Matrix2x3d res;

        Eigen::Vector3d t;
		t = m.tail<3>() + m.head<3>();

        const Eigen::Vector3d T = M - t;
        P = R.transpose() * T;

        const double norm2 = 1.0 / (P.z() * P.z());
        const Eigen::Vector3d k(P.x() * norm2, P.y() * norm2, 1.0 / P.z());

		res(0, 0) = -P.y() * k.x();
		res(1, 0) = -P.tail<2>().dot(k.tail<2>());

		const Eigen::Vector3d tmp1(R.transpose() * T);
		res.col(1) = tmp1.head<2>() * k.z() - k.head<2>() * tmp1.z();

		const Eigen::Vector3d tmp2(R.row(1).transpose() * T.x() - R.row(0).transpose() * T.y());
		
		// commment this line to make it inline everything (else)
		res.col(2) = tmp2.head<2>() * k.z() - tmp2.z() * k.head<2>();
							
		return res;
}

Attachment: eigen_strong_inline.asm.bz2
Description: application/bzip

Attachment: eigen_inline.asm.bz2
Description: application/bzip



Mail converted by MHonArc 2.6.19+ http://listengine.tuxfamily.org/