aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/src
diff options
context:
space:
mode:
Diffstat (limited to 'unsupported/Eigen/src')
-rwxr-xr-x[-rw-r--r--]unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h27
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h1
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h8
-rw-r--r--unsupported/Eigen/src/MatrixFunctions/MatrixPower.h4
-rw-r--r--unsupported/Eigen/src/SparseExtra/RandomSetter.h6
-rw-r--r--unsupported/Eigen/src/Splines/SplineFitting.h6
6 files changed, 30 insertions, 22 deletions
diff --git a/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h b/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
index 8b58b512b..481dfa91a 100644..100755
--- a/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
+++ b/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
@@ -99,7 +99,11 @@ class AutoDiffScalar
{}
template<typename OtherDerType>
- AutoDiffScalar(const AutoDiffScalar<OtherDerType>& other)
+ AutoDiffScalar(const AutoDiffScalar<OtherDerType>& other
+#ifndef EIGEN_PARSED_BY_DOXYGEN
+ , typename internal::enable_if<internal::is_same<Scalar,typename OtherDerType::Scalar>::value,void*>::type = 0
+#endif
+ )
: m_value(other.value()), m_derivatives(other.derivatives())
{}
@@ -127,6 +131,14 @@ class AutoDiffScalar
return *this;
}
+ inline AutoDiffScalar& operator=(const Scalar& other)
+ {
+ m_value = other;
+ if(m_derivatives.size()>0)
+ m_derivatives.setZero();
+ return *this;
+ }
+
// inline operator const Scalar& () const { return m_value; }
// inline operator Scalar& () { return m_value; }
@@ -577,23 +589,24 @@ EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
return ReturnType(log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
template<typename DerType>
-inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<DerType>::Scalar>, const DerType> >
-pow(const Eigen::AutoDiffScalar<DerType>& x, typename Eigen::internal::traits<DerType>::Scalar y)
+inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename internal::traits<typename internal::remove_all<DerType>::type>::Scalar>, const typename internal::remove_all<DerType>::type> >
+pow(const Eigen::AutoDiffScalar<DerType>& x, typename internal::traits<typename internal::remove_all<DerType>::type>::Scalar y)
{
using namespace Eigen;
- typedef typename Eigen::internal::traits<DerType>::Scalar Scalar;
- return AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const DerType> >(
+ typedef typename internal::remove_all<DerType>::type DerTypeCleaned;
+ typedef typename Eigen::internal::traits<DerTypeCleaned>::Scalar Scalar;
+ return AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const DerTypeCleaned> >(
std::pow(x.value(),y),
x.derivatives() * (y * std::pow(x.value(),y-1)));
}
template<typename DerTypeA,typename DerTypeB>
-inline const AutoDiffScalar<Matrix<typename internal::traits<DerTypeA>::Scalar,Dynamic,1> >
+inline const AutoDiffScalar<Matrix<typename internal::traits<typename internal::remove_all<DerTypeA>::type>::Scalar,Dynamic,1> >
atan2(const AutoDiffScalar<DerTypeA>& a, const AutoDiffScalar<DerTypeB>& b)
{
using std::atan2;
- typedef typename internal::traits<DerTypeA>::Scalar Scalar;
+ typedef typename internal::traits<typename internal::remove_all<DerTypeA>::type>::Scalar Scalar;
typedef AutoDiffScalar<Matrix<Scalar,Dynamic,1> > PlainADS;
PlainADS ret;
ret.value() = atan2(a.value(), b.value());
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
index 14a8aef58..bbb7e5776 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixExponential.h
@@ -257,7 +257,6 @@ struct matrix_exp_computeUV<MatrixType, long double>
static void run(const MatrixType& arg, MatrixType& U, MatrixType& V, int& squarings)
{
#if LDBL_MANT_DIG == 53 // double precision
-
matrix_exp_computeUV<MatrixType, double>::run(arg, U, V, squarings);
#else
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h b/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
index 463d7be0c..e43e86e90 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixLogarithm.h
@@ -11,10 +11,6 @@
#ifndef EIGEN_MATRIX_LOGARITHM
#define EIGEN_MATRIX_LOGARITHM
-#ifndef M_PI
-#define M_PI 3.141592653589793238462643383279503L
-#endif
-
namespace Eigen {
namespace internal {
@@ -65,8 +61,8 @@ void matrix_log_compute_2x2(const MatrixType& A, MatrixType& result)
else
{
// computation in previous branch is inaccurate if A(1,1) \approx A(0,0)
- int unwindingNumber = static_cast<int>(ceil((imag(logA11 - logA00) - M_PI) / (2*M_PI)));
- result(0,1) = A(0,1) * (numext::log1p(y/A(0,0)) + Scalar(0,2*M_PI*unwindingNumber)) / y;
+ int unwindingNumber = static_cast<int>(ceil((imag(logA11 - logA00) - EIGEN_PI) / (2*EIGEN_PI)));
+ result(0,1) = A(0,1) * (numext::log1p(y/A(0,0)) + Scalar(0,2*EIGEN_PI*unwindingNumber)) / y;
}
}
diff --git a/unsupported/Eigen/src/MatrixFunctions/MatrixPower.h b/unsupported/Eigen/src/MatrixFunctions/MatrixPower.h
index 1e5a59c55..f37d31c3f 100644
--- a/unsupported/Eigen/src/MatrixFunctions/MatrixPower.h
+++ b/unsupported/Eigen/src/MatrixFunctions/MatrixPower.h
@@ -298,8 +298,8 @@ MatrixPowerAtomic<MatrixType>::computeSuperDiag(const ComplexScalar& curr, const
ComplexScalar logCurr = log(curr);
ComplexScalar logPrev = log(prev);
- int unwindingNumber = ceil((numext::imag(logCurr - logPrev) - M_PI) / (2*M_PI));
- ComplexScalar w = numext::log1p((curr-prev)/prev)/RealScalar(2) + ComplexScalar(0, M_PI*unwindingNumber);
+ int unwindingNumber = ceil((numext::imag(logCurr - logPrev) - EIGEN_PI) / (2*EIGEN_PI));
+ ComplexScalar w = numext::log1p((curr-prev)/prev)/RealScalar(2) + ComplexScalar(0, EIGEN_PI*unwindingNumber);
return RealScalar(2) * exp(RealScalar(0.5) * p * (logCurr + logPrev)) * sinh(p * w) / (curr - prev);
}
diff --git a/unsupported/Eigen/src/SparseExtra/RandomSetter.h b/unsupported/Eigen/src/SparseExtra/RandomSetter.h
index eb3e17330..ee97299af 100644
--- a/unsupported/Eigen/src/SparseExtra/RandomSetter.h
+++ b/unsupported/Eigen/src/SparseExtra/RandomSetter.h
@@ -95,10 +95,10 @@ template<typename Scalar> struct GoogleSparseHashMapTraits
*
* \brief The RandomSetter is a wrapper object allowing to set/update a sparse matrix with random access
*
- * \param SparseMatrixType the type of the sparse matrix we are updating
- * \param MapTraits a traits class representing the map implementation used for the temporary sparse storage.
+ * \tparam SparseMatrixType the type of the sparse matrix we are updating
+ * \tparam MapTraits a traits class representing the map implementation used for the temporary sparse storage.
* Its default value depends on the system.
- * \param OuterPacketBits defines the number of rows (or columns) manage by a single map object
+ * \tparam OuterPacketBits defines the number of rows (or columns) manage by a single map object
* as a power of two exponent.
*
* This class temporarily represents a sparse matrix object using a generic map implementation allowing for
diff --git a/unsupported/Eigen/src/Splines/SplineFitting.h b/unsupported/Eigen/src/Splines/SplineFitting.h
index d3c245fa9..c761a9b3d 100644
--- a/unsupported/Eigen/src/Splines/SplineFitting.h
+++ b/unsupported/Eigen/src/Splines/SplineFitting.h
@@ -130,12 +130,12 @@ namespace Eigen
ParameterVectorType temporaryParameters(numParameters + 1);
KnotVectorType derivativeKnots(numInternalDerivatives);
- for (unsigned int i = 0; i < numAverageKnots - 1; ++i)
+ for (DenseIndex i = 0; i < numAverageKnots - 1; ++i)
{
temporaryParameters[0] = averageKnots[i];
ParameterVectorType parameterIndices(numParameters);
int temporaryParameterIndex = 1;
- for (int j = 0; j < numParameters; ++j)
+ for (DenseIndex j = 0; j < numParameters; ++j)
{
Scalar parameter = parameters[j];
if (parameter >= averageKnots[i] && parameter < averageKnots[i + 1])
@@ -167,7 +167,7 @@ namespace Eigen
derivativeKnots.data(), derivativeKnots.data() + derivativeKnots.size(),
temporaryKnots.data());
- // Number of control points (one for each point and derivative) plus spline order.
+ // Number of knots (one for each point and derivative) plus spline order.
DenseIndex numKnots = numParameters + numDerivatives + degree + 1;
knots.resize(numKnots);