aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-10-16 13:22:38 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-10-16 13:22:38 +0200
commit7b0c4102facc9b5f6ca99ef76febb05a9499b8b0 (patch)
tree074bab8e35f62af40d0e60901edcedf7efb9b411
parent44ba4b1d6d5cd39d824bb83876175d0dc39a9cc3 (diff)
* add a Make* expression type builder to allow the
construction of generic expressions working for both dense and sparse matrix. A nicer solution would be to use CwiseBinaryOp for any kind of matrix. To this end we either need to change the overall design so that the base class(es) depends on the kind of matrix, or we could add a template parameter to each expression type (e.g., int Kind = ei_traits<MatrixType>::Kind) allowing to specialize each expression for each kind of matrix. * Extend AutoDiffScalar to work with sparse vector expression for the derivatives.
-rw-r--r--Eigen/Core1
-rw-r--r--Eigen/Sparse1
-rw-r--r--Eigen/src/Core/ExpressionMaker.h61
-rw-r--r--Eigen/src/Sparse/SparseExpressionMaker.h48
-rw-r--r--unsupported/Eigen/src/AutoDiff/AutoDiffJacobian.h4
-rw-r--r--unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h85
6 files changed, 161 insertions, 39 deletions
diff --git a/Eigen/Core b/Eigen/Core
index c8fcb1c09..3dce6422f 100644
--- a/Eigen/Core
+++ b/Eigen/Core
@@ -200,6 +200,7 @@ namespace Eigen {
#include "src/Core/products/TriangularMatrixMatrix.h"
#include "src/Core/products/TriangularSolverMatrix.h"
#include "src/Core/BandMatrix.h"
+#include "src/Core/ExpressionMaker.h"
} // namespace Eigen
diff --git a/Eigen/Sparse b/Eigen/Sparse
index a8888daa3..96bd61419 100644
--- a/Eigen/Sparse
+++ b/Eigen/Sparse
@@ -110,6 +110,7 @@ namespace Eigen {
#include "src/Sparse/SparseLLT.h"
#include "src/Sparse/SparseLDLT.h"
#include "src/Sparse/SparseLU.h"
+#include "src/Sparse/SparseExpressionMaker.h"
#ifdef EIGEN_CHOLMOD_SUPPORT
# include "src/Sparse/CholmodSupport.h"
diff --git a/Eigen/src/Core/ExpressionMaker.h b/Eigen/src/Core/ExpressionMaker.h
new file mode 100644
index 000000000..1d265b63c
--- /dev/null
+++ b/Eigen/src/Core/ExpressionMaker.h
@@ -0,0 +1,61 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#ifndef EIGEN_EXPRESSIONMAKER_H
+#define EIGEN_EXPRESSIONMAKER_H
+
+// computes the shape of a matrix from its traits flag
+template<typename XprType> struct ei_shape_of
+{
+ enum { ret = ei_traits<XprType>::Flags&SparseBit ? IsSparse : IsDense };
+};
+
+
+// Since the Sparse module is completely separated from the Core module, there is
+// no way to write the type of a generic expression working for both dense and sparse
+// matrix. Unless we change the overall design, here is a workaround.
+// There is an example in unsuported/Eigen/src/AutoDiff/AutoDiffScalar.
+
+template<typename XprType, int Shape = ei_shape_of<XprType>::ret>
+struct MakeNestByValue
+{
+ typedef NestByValue<XprType> Type;
+};
+
+template<typename Func, typename XprType, int Shape = ei_shape_of<XprType>::ret>
+struct MakeCwiseUnaryOp
+{
+ typedef CwiseUnaryOp<Func,XprType> Type;
+};
+
+template<typename Func, typename A, typename B, int Shape = ei_shape_of<A>::ret>
+struct MakeCwiseBinaryOp
+{
+ typedef CwiseBinaryOp<Func,A,B> Type;
+};
+
+// TODO complete the list
+
+
+#endif // EIGEN_EXPRESSIONMAKER_H
diff --git a/Eigen/src/Sparse/SparseExpressionMaker.h b/Eigen/src/Sparse/SparseExpressionMaker.h
new file mode 100644
index 000000000..1fdcbb1f2
--- /dev/null
+++ b/Eigen/src/Sparse/SparseExpressionMaker.h
@@ -0,0 +1,48 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#ifndef EIGEN_SPARSE_EXPRESSIONMAKER_H
+#define EIGEN_SPARSE_EXPRESSIONMAKER_H
+
+template<typename XprType>
+struct MakeNestByValue<XprType,IsSparse>
+{
+ typedef SparseNestByValue<XprType> Type;
+};
+
+template<typename Func, typename XprType>
+struct MakeCwiseUnaryOp<Func,XprType,IsSparse>
+{
+ typedef SparseCwiseUnaryOp<Func,XprType> Type;
+};
+
+template<typename Func, typename A, typename B>
+struct MakeCwiseBinaryOp<Func,A,B,IsSparse>
+{
+ typedef SparseCwiseBinaryOp<Func,A,B> Type;
+};
+
+// TODO complete the list
+
+#endif // EIGEN_SPARSE_EXPRESSIONMAKER_H
diff --git a/unsupported/Eigen/src/AutoDiff/AutoDiffJacobian.h b/unsupported/Eigen/src/AutoDiff/AutoDiffJacobian.h
index a5e881487..b3983f8a6 100644
--- a/unsupported/Eigen/src/AutoDiff/AutoDiffJacobian.h
+++ b/unsupported/Eigen/src/AutoDiff/AutoDiffJacobian.h
@@ -50,10 +50,12 @@ public:
typedef typename Functor::InputType InputType;
typedef typename Functor::ValueType ValueType;
typedef typename Functor::JacobianType JacobianType;
+ typedef typename JacobianType::Scalar Scalar;
- typedef Matrix<double,InputsAtCompileTime,1> DerivativeType;
+ typedef Matrix<Scalar,InputsAtCompileTime,1> DerivativeType;
typedef AutoDiffScalar<DerivativeType> ActiveScalar;
+
typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
diff --git a/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h b/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
index fc5e237ab..2fb733a99 100644
--- a/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
+++ b/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
@@ -42,9 +42,17 @@ void ei_make_coherent(const A& a, const B&b)
/** \class AutoDiffScalar
* \brief A scalar type replacement with automatic differentation capability
*
- * \param DerType the vector type used to store/represent the derivatives (e.g. Vector3f)
+ * \param _DerType the vector type used to store/represent the derivatives. The base scalar type
+ * as well as the number of derivatives to compute are determined from this type.
+ * Typical choices include, e.g., \c Vector4f for 4 derivatives, or \c VectorXf
+ * if the number of derivatives is not known at compile time, and/or, the number
+ * of derivatives is large.
+ * Note that _DerType can also be a reference (e.g., \c VectorXf&) to wrap a
+ * existing vector into an AutoDiffScalar.
+ * Finally, _DerType can also be any Eigen compatible expression.
*
- * This class represents a scalar value while tracking its respective derivatives.
+ * This class represents a scalar value while tracking its respective derivatives using Eigen's expression
+ * template mechanism.
*
* It supports the following list of global math function:
* - std::abs, std::sqrt, std::pow, std::exp, std::log, std::sin, std::cos,
@@ -56,10 +64,11 @@ void ei_make_coherent(const A& a, const B&b)
* while derivatives are computed right away.
*
*/
-template<typename DerType>
+template<typename _DerType>
class AutoDiffScalar
{
public:
+ typedef typename ei_cleantype<_DerType>::type DerType;
typedef typename ei_traits<DerType>::Scalar Scalar;
inline AutoDiffScalar() {}
@@ -108,12 +117,12 @@ class AutoDiffScalar
inline const DerType& derivatives() const { return m_derivatives; }
inline DerType& derivatives() { return m_derivatives; }
- inline const AutoDiffScalar<DerType> operator+(const Scalar& other) const
+ inline const AutoDiffScalar<DerType&> operator+(const Scalar& other) const
{
return AutoDiffScalar<DerType>(m_value + other, m_derivatives);
}
- friend inline const AutoDiffScalar<DerType> operator+(const Scalar& a, const AutoDiffScalar& b)
+ friend inline const AutoDiffScalar<DerType&> operator+(const Scalar& a, const AutoDiffScalar& b)
{
return AutoDiffScalar<DerType>(a + b.value(), b.derivatives());
}
@@ -125,11 +134,11 @@ class AutoDiffScalar
}
template<typename OtherDerType>
- inline const AutoDiffScalar<CwiseBinaryOp<ei_scalar_sum_op<Scalar>,DerType,OtherDerType> >
+ inline const AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_sum_op<Scalar>,DerType,typename ei_cleantype<OtherDerType>::type>::Type >
operator+(const AutoDiffScalar<OtherDerType>& other) const
{
ei_make_coherent(m_derivatives, other.derivatives());
- return AutoDiffScalar<CwiseBinaryOp<ei_scalar_sum_op<Scalar>,DerType,OtherDerType> >(
+ return AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_sum_op<Scalar>,DerType,typename ei_cleantype<OtherDerType>::type>::Type >(
m_value + other.value(),
m_derivatives + other.derivatives());
}
@@ -143,11 +152,11 @@ class AutoDiffScalar
}
template<typename OtherDerType>
- inline const AutoDiffScalar<CwiseBinaryOp<ei_scalar_difference_op<Scalar>, DerType,OtherDerType> >
+ inline const AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_difference_op<Scalar>, DerType,typename ei_cleantype<OtherDerType>::type>::Type >
operator-(const AutoDiffScalar<OtherDerType>& other) const
{
ei_make_coherent(m_derivatives, other.derivatives());
- return AutoDiffScalar<CwiseBinaryOp<ei_scalar_difference_op<Scalar>, DerType,OtherDerType> >(
+ return AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_difference_op<Scalar>, DerType,typename ei_cleantype<OtherDerType>::type>::Type >(
m_value - other.value(),
m_derivatives - other.derivatives());
}
@@ -161,73 +170,73 @@ class AutoDiffScalar
}
template<typename OtherDerType>
- inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_opposite_op<Scalar>, DerType> >
+ inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_opposite_op<Scalar>, DerType>::Type >
operator-() const
{
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_opposite_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_opposite_op<Scalar>, DerType>::Type >(
-m_value,
-m_derivatives);
}
- inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >
+ inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >
operator*(const Scalar& other) const
{
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >(
m_value * other,
(m_derivatives * other));
}
- friend inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >
+ friend inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >
operator*(const Scalar& other, const AutoDiffScalar& a)
{
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >(
a.value() * other,
a.derivatives() * other);
}
- inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >
+ inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >
operator/(const Scalar& other) const
{
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >(
m_value / other,
(m_derivatives * (Scalar(1)/other)));
}
- friend inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >
+ friend inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >
operator/(const Scalar& other, const AutoDiffScalar& a)
{
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >(
other / a.value(),
a.derivatives() * (-Scalar(1)/other));
}
template<typename OtherDerType>
- inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>,
- NestByValue<CwiseBinaryOp<ei_scalar_difference_op<Scalar>,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, OtherDerType> > > > > >
+ inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseBinaryOp<ei_scalar_difference_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type>::Type,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, typename ei_cleantype<OtherDerType>::type>::Type>::Type >::Type >::Type >::Type >
operator/(const AutoDiffScalar<OtherDerType>& other) const
{
ei_make_coherent(m_derivatives, other.derivatives());
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>,
- NestByValue<CwiseBinaryOp<ei_scalar_difference_op<Scalar>,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, OtherDerType> > > > > >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseBinaryOp<ei_scalar_difference_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type>::Type,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, typename ei_cleantype<OtherDerType>::type>::Type>::Type >::Type >::Type >::Type >(
m_value / other.value(),
((m_derivatives * other.value()).nestByValue() - (m_value * other.derivatives()).nestByValue()).nestByValue()
* (Scalar(1)/(other.value()*other.value())));
}
template<typename OtherDerType>
- inline const AutoDiffScalar<CwiseBinaryOp<ei_scalar_sum_op<Scalar>,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, OtherDerType> > > >
+ inline const AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_sum_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type>::Type,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, typename ei_cleantype<OtherDerType>::type>::Type>::Type >::Type >
operator*(const AutoDiffScalar<OtherDerType>& other) const
{
ei_make_coherent(m_derivatives, other.derivatives());
- return AutoDiffScalar<CwiseBinaryOp<ei_scalar_sum_op<Scalar>,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >,
- NestByValue<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, OtherDerType> > > >(
+ return AutoDiffScalar<typename MakeCwiseBinaryOp<ei_scalar_sum_op<Scalar>,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type>::Type,
+ typename MakeNestByValue<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, typename ei_cleantype<OtherDerType>::type>::Type>::Type >::Type >(
m_value * other.value(),
(m_derivatives * other.value()).nestByValue() + (m_value * other.derivatives()).nestByValue());
}
@@ -299,11 +308,11 @@ struct ei_make_coherent_impl<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRo
#define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
template<typename DerType> \
- inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::ei_scalar_multiple_op<typename Eigen::ei_traits<DerType>::Scalar>, DerType> > \
+ inline const Eigen::AutoDiffScalar<typename Eigen::MakeCwiseUnaryOp<Eigen::ei_scalar_multiple_op<typename Eigen::ei_traits<DerType>::Scalar>, DerType>::Type > \
FUNC(const Eigen::AutoDiffScalar<DerType>& x) { \
using namespace Eigen; \
typedef typename ei_traits<DerType>::Scalar Scalar; \
- typedef AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> > ReturnType; \
+ typedef AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type > ReturnType; \
CODE; \
}
@@ -330,12 +339,12 @@ namespace std
return ReturnType(std::log(x.value),x.derivatives() * (Scalar(1).x.value()));)
template<typename DerType>
- inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::ei_scalar_multiple_op<typename Eigen::ei_traits<DerType>::Scalar>, DerType> >
+ inline const Eigen::AutoDiffScalar<typename Eigen::MakeCwiseUnaryOp<Eigen::ei_scalar_multiple_op<typename Eigen::ei_traits<DerType>::Scalar>, DerType>::Type >
pow(const Eigen::AutoDiffScalar<DerType>& x, typename Eigen::ei_traits<DerType>::Scalar y)
{
using namespace Eigen;
typedef typename ei_traits<DerType>::Scalar Scalar;
- return AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType> >(
+ return AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<Scalar>, DerType>::Type >(
std::pow(x.value(),y),
x.derivatives() * (y * std::pow(x.value(),y-1)));
}
@@ -375,7 +384,7 @@ EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(ei_log,
return ReturnType(ei_log(x.value),x.derivatives() * (Scalar(1).x.value()));)
template<typename DerType>
-inline const AutoDiffScalar<CwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<DerType>::Scalar>, DerType> >
+inline const AutoDiffScalar<typename MakeCwiseUnaryOp<ei_scalar_multiple_op<typename ei_traits<DerType>::Scalar>, DerType>::Type >
ei_pow(const AutoDiffScalar<DerType>& x, typename ei_traits<DerType>::Scalar y)
{ return std::pow(x,y);}