aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Core/Assign.h2
-rw-r--r--Eigen/src/Core/Block.h26
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h13
-rw-r--r--Eigen/src/Core/CwiseNullaryOp.h4
-rw-r--r--Eigen/src/Core/Dot.h6
-rw-r--r--Eigen/src/Core/Functors.h2
-rw-r--r--Eigen/src/Core/Fuzzy.h8
-rw-r--r--Eigen/src/Core/Map.h2
-rw-r--r--Eigen/src/Core/Matrix.h8
-rw-r--r--Eigen/src/Core/Part.h4
-rw-r--r--Eigen/src/Core/Product.h14
-rw-r--r--Eigen/src/Core/util/Constants.h4
-rw-r--r--Eigen/src/Core/util/StaticAssert.h15
-rw-r--r--Eigen/src/Geometry/Hyperplane.h4
-rw-r--r--Eigen/src/Geometry/OrthoMethods.h4
-rw-r--r--Eigen/src/Geometry/ParametrizedLine.h2
-rw-r--r--Eigen/src/Geometry/Rotation2D.h2
-rw-r--r--Eigen/src/Geometry/RotationBase.h8
-rw-r--r--Eigen/src/Geometry/Transform.h20
-rw-r--r--Eigen/src/LU/Inverse.h2
-rw-r--r--Eigen/src/Sparse/SparseBlock.h2
-rw-r--r--test/regression.cpp6
22 files changed, 88 insertions, 70 deletions
diff --git a/Eigen/src/Core/Assign.h b/Eigen/src/Core/Assign.h
index fb4857803..fffa78919 100644
--- a/Eigen/src/Core/Assign.h
+++ b/Eigen/src/Core/Assign.h
@@ -400,7 +400,7 @@ template<typename OtherDerived>
inline Derived& MatrixBase<Derived>
::lazyAssign(const MatrixBase<OtherDerived>& other)
{
- EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived);
+ EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived)
ei_assert(rows() == other.rows() && cols() == other.cols());
ei_assign_impl<Derived, OtherDerived>::run(derived(),other.derived());
return derived();
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index 8661352dd..b9a727080 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -122,7 +122,7 @@ template<typename MatrixType, int BlockRows, int BlockCols, int PacketAccess, in
: m_matrix(matrix), m_startRow(startRow), m_startCol(startCol),
m_blockRows(matrix.rows()), m_blockCols(matrix.cols())
{
- EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic,this_method_is_only_for_fixed_size);
+ EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic,this_method_is_only_for_fixed_size)
ei_assert(startRow >= 0 && BlockRows >= 1 && startRow + BlockRows <= matrix.rows()
&& startCol >= 0 && BlockCols >= 1 && startCol + BlockCols <= matrix.cols());
}
@@ -340,7 +340,7 @@ template<typename Derived>
inline typename BlockReturnType<Derived>::SubVectorType MatrixBase<Derived>
::segment(int start, int size)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return typename BlockReturnType<Derived>::SubVectorType(derived(), RowsAtCompileTime == 1 ? 0 : start,
ColsAtCompileTime == 1 ? 0 : start,
RowsAtCompileTime == 1 ? 1 : size,
@@ -352,7 +352,7 @@ template<typename Derived>
inline const typename BlockReturnType<Derived>::SubVectorType
MatrixBase<Derived>::segment(int start, int size) const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return typename BlockReturnType<Derived>::SubVectorType(derived(), RowsAtCompileTime == 1 ? 0 : start,
ColsAtCompileTime == 1 ? 0 : start,
RowsAtCompileTime == 1 ? 1 : size,
@@ -380,7 +380,7 @@ template<typename Derived>
inline typename BlockReturnType<Derived,Dynamic>::SubVectorType
MatrixBase<Derived>::start(int size)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived,
RowsAtCompileTime == 1 ? 1 : Dynamic,
ColsAtCompileTime == 1 ? 1 : Dynamic>
@@ -394,7 +394,7 @@ template<typename Derived>
inline const typename BlockReturnType<Derived,Dynamic>::SubVectorType
MatrixBase<Derived>::start(int size) const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived,
RowsAtCompileTime == 1 ? 1 : Dynamic,
ColsAtCompileTime == 1 ? 1 : Dynamic>
@@ -424,7 +424,7 @@ template<typename Derived>
inline typename BlockReturnType<Derived,Dynamic>::SubVectorType
MatrixBase<Derived>::end(int size)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived,
RowsAtCompileTime == 1 ? 1 : Dynamic,
ColsAtCompileTime == 1 ? 1 : Dynamic>
@@ -440,7 +440,7 @@ template<typename Derived>
inline const typename BlockReturnType<Derived,Dynamic>::SubVectorType
MatrixBase<Derived>::end(int size) const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived,
RowsAtCompileTime == 1 ? 1 : Dynamic,
ColsAtCompileTime == 1 ? 1 : Dynamic>
@@ -469,7 +469,7 @@ template<int Size>
inline typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::segment(int start)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, (RowsAtCompileTime == 1 ? 1 : Size),
(ColsAtCompileTime == 1 ? 1 : Size)>
(derived(), RowsAtCompileTime == 1 ? 0 : start,
@@ -482,7 +482,7 @@ template<int Size>
inline const typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::segment(int start) const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, (RowsAtCompileTime == 1 ? 1 : Size),
(ColsAtCompileTime == 1 ? 1 : Size)>
(derived(), RowsAtCompileTime == 1 ? 0 : start,
@@ -507,7 +507,7 @@ template<int Size>
inline typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::start()
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, (RowsAtCompileTime == 1 ? 1 : Size),
(ColsAtCompileTime == 1 ? 1 : Size)>(derived(), 0, 0);
}
@@ -518,7 +518,7 @@ template<int Size>
inline const typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::start() const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, (RowsAtCompileTime == 1 ? 1 : Size),
(ColsAtCompileTime == 1 ? 1 : Size)>(derived(), 0, 0);
}
@@ -539,7 +539,7 @@ template<int Size>
inline typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::end()
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, RowsAtCompileTime == 1 ? 1 : Size,
ColsAtCompileTime == 1 ? 1 : Size>
(derived(),
@@ -553,7 +553,7 @@ template<int Size>
inline const typename BlockReturnType<Derived,Size>::SubVectorType
MatrixBase<Derived>::end() const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return Block<Derived, RowsAtCompileTime == 1 ? 1 : Size,
ColsAtCompileTime == 1 ? 1 : Size>
(derived(),
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index e28d4266f..c2dad094f 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -46,12 +46,15 @@
template<typename BinaryOp, typename Lhs, typename Rhs>
struct ei_traits<CwiseBinaryOp<BinaryOp, Lhs, Rhs> >
{
+ // even though we require Lhs and Rhs to have the same scalar type (see CwiseBinaryOp constructor),
+ // we still want to handle the case when the result type is different.
typedef typename ei_result_of<
BinaryOp(
typename Lhs::Scalar,
typename Rhs::Scalar
)
>::type Scalar;
+
typedef typename Lhs::Nested LhsNested;
typedef typename Rhs::Nested RhsNested;
typedef typename ei_unref<LhsNested>::type _LhsNested;
@@ -89,6 +92,16 @@ class CwiseBinaryOp : ei_no_assignment_operator,
inline CwiseBinaryOp(const Lhs& lhs, const Rhs& rhs, const BinaryOp& func = BinaryOp())
: m_lhs(lhs), m_rhs(rhs), m_functor(func)
{
+ // we require Lhs and Rhs to have the same scalar type. Currently there is no example of a binary functor
+ // that would take two operands of different types. If there were such an example, then this check should then be
+ // moved to the BinaryOp functors, on a per-case basis. This would however require a change in the BinaryOp functors, as
+ // currently they take only one typename Scalar template parameter.
+ // It is tempting to always allow mixing different types but remember that this is often impossible in the vectorized paths.
+ // So allowing mixing different types gives very unexpected errors when enabling vectorization, when the user tries to
+ // add together a float matrix and a double matrix.
+ EIGEN_STATIC_ASSERT((ei_is_same_type<typename Lhs::Scalar, typename Rhs::Scalar>::ret),
+ you_mixed_different_numeric_types__you_need_to_use_the_cast_method_of_MatrixBase_to_cast_numeric_types_explicitly)
+ // require the sizes to be the same
ei_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
}
diff --git a/Eigen/src/Core/CwiseNullaryOp.h b/Eigen/src/Core/CwiseNullaryOp.h
index 005340213..fb63ad4fd 100644
--- a/Eigen/src/Core/CwiseNullaryOp.h
+++ b/Eigen/src/Core/CwiseNullaryOp.h
@@ -561,7 +561,7 @@ inline Derived& MatrixBase<Derived>::setIdentity()
template<typename Derived>
const typename MatrixBase<Derived>::BasisReturnType MatrixBase<Derived>::Unit(int size, int i)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return BasisReturnType(SquareMatrixType::Identity(size,size), i);
}
@@ -576,7 +576,7 @@ const typename MatrixBase<Derived>::BasisReturnType MatrixBase<Derived>::Unit(in
template<typename Derived>
const typename MatrixBase<Derived>::BasisReturnType MatrixBase<Derived>::Unit(int i)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return BasisReturnType(SquareMatrixType::Identity(),i);
}
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index 3032f79ec..e700b76ae 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -263,9 +263,9 @@ MatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
typedef typename ei_unref<Nested>::type _Nested;
typedef typename ei_unref<OtherNested>::type _OtherNested;
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(_Nested);
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(_OtherNested);
- EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(_Nested,_OtherNested);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(_Nested)
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(_OtherNested)
+ EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(_Nested,_OtherNested)
ei_assert(size() == other.size());
return ei_dot_impl<_Nested, _OtherNested>::run(derived(), other.derived());
diff --git a/Eigen/src/Core/Functors.h b/Eigen/src/Core/Functors.h
index d68558e7b..de45c294d 100644
--- a/Eigen/src/Core/Functors.h
+++ b/Eigen/src/Core/Functors.h
@@ -348,7 +348,7 @@ template<typename Scalar>
struct ei_functor_traits<ei_scalar_identity_op<Scalar> >
{ enum { Cost = NumTraits<Scalar>::AddCost, PacketAccess = false, IsRepeatable = true }; };
-// NOTE quick hack:
+// FIXME quick hack:
// all functors allow linear access, except ei_scalar_identity_op. So we fix here a quick meta
// to indicate whether a functor allows linear access, just always answering 'yes' except for
// ei_scalar_identity_op.
diff --git a/Eigen/src/Core/Fuzzy.h b/Eigen/src/Core/Fuzzy.h
index 7c2acb688..18150cc6d 100644
--- a/Eigen/src/Core/Fuzzy.h
+++ b/Eigen/src/Core/Fuzzy.h
@@ -176,7 +176,7 @@ struct ei_fuzzy_selector<Derived,OtherDerived,true>
typedef typename Derived::RealScalar RealScalar;
static bool isApprox(const Derived& self, const OtherDerived& other, RealScalar prec)
{
- EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived);
+ EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
ei_assert(self.size() == other.size());
return((self - other).squaredNorm() <= std::min(self.squaredNorm(), other.squaredNorm()) * prec * prec);
}
@@ -186,7 +186,7 @@ struct ei_fuzzy_selector<Derived,OtherDerived,true>
}
static bool isMuchSmallerThan(const Derived& self, const OtherDerived& other, RealScalar prec)
{
- EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived);
+ EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
ei_assert(self.size() == other.size());
return(self.squaredNorm() <= other.squaredNorm() * prec * prec);
}
@@ -198,7 +198,7 @@ struct ei_fuzzy_selector<Derived,OtherDerived,false>
typedef typename Derived::RealScalar RealScalar;
static bool isApprox(const Derived& self, const OtherDerived& other, RealScalar prec)
{
- EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived);
+ EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived)
ei_assert(self.rows() == other.rows() && self.cols() == other.cols());
typename Derived::Nested nested(self);
typename OtherDerived::Nested otherNested(other);
@@ -218,7 +218,7 @@ struct ei_fuzzy_selector<Derived,OtherDerived,false>
}
static bool isMuchSmallerThan(const Derived& self, const OtherDerived& other, RealScalar prec)
{
- EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived);
+ EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,OtherDerived)
ei_assert(self.rows() == other.rows() && self.cols() == other.cols());
typename Derived::Nested nested(self);
typename OtherDerived::Nested otherNested(other);
diff --git a/Eigen/src/Core/Map.h b/Eigen/src/Core/Map.h
index fbee1214c..d391ba93f 100644
--- a/Eigen/src/Core/Map.h
+++ b/Eigen/src/Core/Map.h
@@ -90,7 +90,7 @@ template<typename MatrixType, int PacketAccess> class Map
inline void resize(int size)
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(MatrixType);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(MatrixType)
EIGEN_ONLY_USED_FOR_DEBUG(size);
ei_assert(size == this->size());
}
diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h
index 86205aab2..8981f0853 100644
--- a/Eigen/src/Core/Matrix.h
+++ b/Eigen/src/Core/Matrix.h
@@ -369,21 +369,21 @@ class Matrix
/** constructs an initialized 2D vector with given coefficients */
inline Matrix(const float& x, const float& y)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 2);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 2)
m_storage.data()[0] = x;
m_storage.data()[1] = y;
}
/** constructs an initialized 2D vector with given coefficients */
inline Matrix(const double& x, const double& y)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 2);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 2)
m_storage.data()[0] = x;
m_storage.data()[1] = y;
}
/** constructs an initialized 3D vector with given coefficients */
inline Matrix(const Scalar& x, const Scalar& y, const Scalar& z)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 3);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 3)
m_storage.data()[0] = x;
m_storage.data()[1] = y;
m_storage.data()[2] = z;
@@ -391,7 +391,7 @@ class Matrix
/** constructs an initialized 4D vector with given coefficients */
inline Matrix(const Scalar& x, const Scalar& y, const Scalar& z, const Scalar& w)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 4);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Matrix, 4)
m_storage.data()[0] = x;
m_storage.data()[1] = y;
m_storage.data()[2] = z;
diff --git a/Eigen/src/Core/Part.h b/Eigen/src/Core/Part.h
index b7daf4e4e..22eb65318 100644
--- a/Eigen/src/Core/Part.h
+++ b/Eigen/src/Core/Part.h
@@ -100,8 +100,8 @@ template<typename MatrixType, unsigned int Mode> class Part
inline Scalar& coeffRef(int row, int col)
{
- EIGEN_STATIC_ASSERT(!(Flags & UnitDiagBit), writting_to_triangular_part_with_unit_diag_is_not_supported);
- EIGEN_STATIC_ASSERT(!(Flags & SelfAdjointBit), default_writting_to_selfadjoint_not_supported);
+ EIGEN_STATIC_ASSERT(!(Flags & UnitDiagBit), writing_to_triangular_part_with_unit_diagonal_is_not_supported)
+ EIGEN_STATIC_ASSERT(!(Flags & SelfAdjointBit), default_writing_to_selfadjoint_not_supported)
ei_assert( (Mode==Upper && col>=row)
|| (Mode==Lower && col<=row)
|| (Mode==StrictlyUpper && col>row)
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index ea9d739fa..c5b06c450 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -116,8 +116,8 @@ template<typename LhsNested, typename RhsNested, int ProductMode>
struct ei_traits<Product<LhsNested, RhsNested, ProductMode> >
{
// clean the nested types:
- typedef typename ei_unconst<typename ei_unref<LhsNested>::type>::type _LhsNested;
- typedef typename ei_unconst<typename ei_unref<RhsNested>::type>::type _RhsNested;
+ typedef typename ei_cleantype<LhsNested>::type _LhsNested;
+ typedef typename ei_cleantype<RhsNested>::type _RhsNested;
typedef typename _LhsNested::Scalar Scalar;
enum {
@@ -194,6 +194,10 @@ template<typename LhsNested, typename RhsNested, int ProductMode> class Product
inline Product(const Lhs& lhs, const Rhs& rhs)
: m_lhs(lhs), m_rhs(rhs)
{
+ // we don't allow taking products of matrices of different real types, as that wouldn't be vectorizable.
+ // We still allow to mix T and complex<T>.
+ EIGEN_STATIC_ASSERT((ei_is_same_type<typename Lhs::RealScalar, typename Rhs::RealScalar>::ret),
+ you_mixed_different_numeric_types__you_need_to_use_the_cast_method_of_MatrixBase_to_cast_numeric_types_explicitly)
ei_assert(lhs.cols() == rhs.rows()
&& "invalid matrix product"
&& "if you wanted a coeff-wise or a dot product use the respective explicit functions");
@@ -278,10 +282,10 @@ MatrixBase<Derived>::operator*(const MatrixBase<OtherDerived> &other) const
// * for a dot product use: v1.dot(v2)
// * for a coeff-wise product use: v1.cwise()*v2
EIGEN_STATIC_ASSERT(ProductIsValid || !(AreVectors && SameSizes),
- invalid_vector_vector_product__if_you_wanted_a_dot_or_coeff_wise_product_you_must_use_the_explicit_functions);
+ invalid_vector_vector_product__if_you_wanted_a_dot_or_coeff_wise_product_you_must_use_the_explicit_functions)
EIGEN_STATIC_ASSERT(ProductIsValid || !(SameSizes && !AreVectors),
- invalid_matrix_product__if_you_wanted_a_coeff_wise_product_you_must_use_the_explicit_function);
- EIGEN_STATIC_ASSERT(ProductIsValid || SameSizes, invalid_matrix_product);
+ invalid_matrix_product__if_you_wanted_a_coeff_wise_product_you_must_use_the_explicit_function)
+ EIGEN_STATIC_ASSERT(ProductIsValid || SameSizes, invalid_matrix_product)
return typename ProductReturnType<Derived,OtherDerived>::Type(derived(), other.derived());
}
diff --git a/Eigen/src/Core/util/Constants.h b/Eigen/src/Core/util/Constants.h
index e83c9f74c..20b427df1 100644
--- a/Eigen/src/Core/util/Constants.h
+++ b/Eigen/src/Core/util/Constants.h
@@ -31,14 +31,14 @@
*
* Explanation for the choice of this value:
* - It should be positive and larger than any reasonable compile-time-fixed number of rows or columns.
- * This means that it should be at least 128 or so.
+ * This allows to simplify many compile-time conditions throughout Eigen.
* - It should be smaller than the sqrt of INT_MAX. Indeed, we often multiply a number of rows with a number
* of columns in order to compute a number of coefficients. Even if we guard that with an "if" checking whether
* the values are Dynamic, we still get a compiler warning "integer overflow". So the only way to get around
* it would be a meta-selector. Doing this everywhere would reduce code readability and lenghten compilation times.
* Also, disabling compiler warnings for integer overflow, sounds like a bad idea.
*
- * If you wish to port Eigen to a platform where sizeof(int)==2, it is perfectly possible to set Dynamic to, say, 250.
+ * If you wish to port Eigen to a platform where sizeof(int)==2, it is perfectly possible to set Dynamic to, say, 100.
*/
const int Dynamic = 10000;
diff --git a/Eigen/src/Core/util/StaticAssert.h b/Eigen/src/Core/util/StaticAssert.h
index 4f22cf5f2..44f0008be 100644
--- a/Eigen/src/Core/util/StaticAssert.h
+++ b/Eigen/src/Core/util/StaticAssert.h
@@ -60,23 +60,24 @@
you_mixed_matrices_of_different_sizes,
this_method_is_only_for_vectors_of_a_specific_size,
this_method_is_only_for_matrices_of_a_specific_size,
- you_did_a_programming_error,
+ you_made_a_programming_mistake,
you_called_a_fixed_size_method_on_a_dynamic_size_matrix_or_vector,
unaligned_load_and_store_operations_unimplemented_on_AltiVec,
- scalar_type_must_be_floating_point,
- default_writting_to_selfadjoint_not_supported,
- writting_to_triangular_part_with_unit_diag_is_not_supported,
+ numeric_type_must_be_floating_point,
+ default_writing_to_selfadjoint_not_supported,
+ writing_to_triangular_part_with_unit_diagonal_is_not_supported,
this_method_is_only_for_fixed_size,
invalid_matrix_product,
invalid_vector_vector_product__if_you_wanted_a_dot_or_coeff_wise_product_you_must_use_the_explicit_functions,
- invalid_matrix_product__if_you_wanted_a_coeff_wise_product_you_must_use_the_explicit_function
+ invalid_matrix_product__if_you_wanted_a_coeff_wise_product_you_must_use_the_explicit_function,
+ you_mixed_different_numeric_types__you_need_to_use_the_cast_method_of_MatrixBase_to_cast_numeric_types_explicitly
};
};
#define EIGEN_STATIC_ASSERT(CONDITION,MSG) \
if (Eigen::ei_static_assert<CONDITION ? true : false>::MSG) {}
- #endif // CXX0X
+ #endif // not CXX0X
#else // EIGEN_NO_STATIC_ASSERT
@@ -121,7 +122,7 @@
|| int(TYPE1::ColsAtCompileTime)==Eigen::Dynamic \
|| int(TYPE0::ColsAtCompileTime)==int(TYPE1::ColsAtCompileTime)))
-// static assertion failing if the two matrix expression types are not compatible (same fixed-size or dynamic size)
+// static assertion failing if it is guaranteed at compile-time that the two matrix expression types have different sizes
#define EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(TYPE0,TYPE1) \
EIGEN_STATIC_ASSERT( \
EIGEN_PREDICATE_SAME_MATRIX_SIZE(TYPE0,TYPE1),\
diff --git a/Eigen/src/Geometry/Hyperplane.h b/Eigen/src/Geometry/Hyperplane.h
index 599579726..5171f9ae4 100644
--- a/Eigen/src/Geometry/Hyperplane.h
+++ b/Eigen/src/Geometry/Hyperplane.h
@@ -104,7 +104,7 @@ public:
*/
static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
Hyperplane result(p0.size());
result.normal() = (p2 - p0).cross(p1 - p0).normalized();
result.offset() = -result.normal().dot(p0);
@@ -184,7 +184,7 @@ public:
*/
VectorType intersection(const Hyperplane& other)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
// since the line equations ax+by=c are normalized with a^2+b^2=1, the following tests
// whether the two lines are approximately parallel.
diff --git a/Eigen/src/Geometry/OrthoMethods.h b/Eigen/src/Geometry/OrthoMethods.h
index f545f176d..d60664cc0 100644
--- a/Eigen/src/Geometry/OrthoMethods.h
+++ b/Eigen/src/Geometry/OrthoMethods.h
@@ -37,7 +37,7 @@ template<typename OtherDerived>
inline typename MatrixBase<Derived>::EvalType
MatrixBase<Derived>::cross(const MatrixBase<OtherDerived>& other) const
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3)
// Note that there is no need for an expression here since the compiler
// optimize such a small temporary very well (even within a complex expression)
@@ -112,7 +112,7 @@ template<typename Derived>
typename MatrixBase<Derived>::EvalType
MatrixBase<Derived>::unitOrthogonal() const
{
- EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
+ EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
return ei_unitOrthogonal_selector<Derived>::run(derived());
}
diff --git a/Eigen/src/Geometry/ParametrizedLine.h b/Eigen/src/Geometry/ParametrizedLine.h
index 761dea1e1..b344493c7 100644
--- a/Eigen/src/Geometry/ParametrizedLine.h
+++ b/Eigen/src/Geometry/ParametrizedLine.h
@@ -141,7 +141,7 @@ protected:
template <typename _Scalar, int _AmbientDim>
inline ParametrizedLine<_Scalar, _AmbientDim>::ParametrizedLine(const Hyperplane<_Scalar, _AmbientDim>& hyperplane)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2);
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
direction() = hyperplane.normal().unitOrthogonal();
origin() = -hyperplane.normal()*hyperplane.offset();
}
diff --git a/Eigen/src/Geometry/Rotation2D.h b/Eigen/src/Geometry/Rotation2D.h
index aac7648b1..6bc7328ea 100644
--- a/Eigen/src/Geometry/Rotation2D.h
+++ b/Eigen/src/Geometry/Rotation2D.h
@@ -140,7 +140,7 @@ template<typename Scalar>
template<typename Derived>
Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
{
- EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,you_made_a_programming_mistake)
m_angle = ei_atan2(mat.coeff(1,0), mat.coeff(0,0));
return *this;
}
diff --git a/Eigen/src/Geometry/RotationBase.h b/Eigen/src/Geometry/RotationBase.h
index dff905528..6ad057e03 100644
--- a/Eigen/src/Geometry/RotationBase.h
+++ b/Eigen/src/Geometry/RotationBase.h
@@ -77,7 +77,7 @@ template<typename OtherDerived>
Matrix<_Scalar, _Rows, _Cols, _Storage, _MaxRows, _MaxCols>
::Matrix(const RotationBase<OtherDerived,ColsAtCompileTime>& r)
{
- EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix,int(OtherDerived::Dim),int(OtherDerived::Dim));
+ EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix,int(OtherDerived::Dim),int(OtherDerived::Dim))
*this = r.toRotationMatrix();
}
@@ -91,7 +91,7 @@ Matrix<_Scalar, _Rows, _Cols, _Storage, _MaxRows, _MaxCols>&
Matrix<_Scalar, _Rows, _Cols, _Storage, _MaxRows, _MaxCols>
::operator=(const RotationBase<OtherDerived,ColsAtCompileTime>& r)
{
- EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix,int(OtherDerived::Dim),int(OtherDerived::Dim));
+ EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix,int(OtherDerived::Dim),int(OtherDerived::Dim))
return *this = r.toRotationMatrix();
}
@@ -116,7 +116,7 @@ Matrix<_Scalar, _Rows, _Cols, _Storage, _MaxRows, _MaxCols>
template<typename Scalar, int Dim>
inline static Matrix<Scalar,2,2> ei_toRotationMatrix(const Scalar& s)
{
- EIGEN_STATIC_ASSERT(Dim==2,you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Dim==2,you_made_a_programming_mistake)
return Rotation2D<Scalar>(s).toRotationMatrix();
}
@@ -130,7 +130,7 @@ template<typename Scalar, int Dim, typename OtherDerived>
inline static const MatrixBase<OtherDerived>& ei_toRotationMatrix(const MatrixBase<OtherDerived>& mat)
{
EIGEN_STATIC_ASSERT(OtherDerived::RowsAtCompileTime==Dim && OtherDerived::ColsAtCompileTime==Dim,
- you_did_a_programming_error);
+ you_made_a_programming_mistake)
return mat;
}
diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h
index a0c07e405..00cc8e47a 100644
--- a/Eigen/src/Geometry/Transform.h
+++ b/Eigen/src/Geometry/Transform.h
@@ -302,7 +302,7 @@ Transform<Scalar,Dim>::Transform(const QMatrix& other)
template<typename Scalar, int Dim>
Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QMatrix& other)
{
- EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Dim==2, you_made_a_programming_mistake)
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy(),
0, 0, 1;
@@ -318,7 +318,7 @@ Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QMatrix& other)
template<typename Scalar, int Dim>
QMatrix Transform<Scalar,Dim>::toQMatrix(void) const
{
- EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Dim==2, you_made_a_programming_mistake)
return QMatrix(other.coeffRef(0,0), other.coeffRef(1,0),
other.coeffRef(0,1), other.coeffRef(1,1),
other.coeffRef(0,2), other.coeffRef(1,2));
@@ -341,7 +341,7 @@ Transform<Scalar,Dim>::Transform(const QTransform& other)
template<typename Scalar, int Dim>
Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QTransform& other)
{
- EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Dim==2, you_made_a_programming_mistake)
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy(),
other.m13(), other.m23(), other.m33();
@@ -355,7 +355,7 @@ Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QTransform& other)
template<typename Scalar, int Dim>
QMatrix Transform<Scalar,Dim>::toQTransform(void) const
{
- EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(Dim==2, you_made_a_programming_mistake)
return QTransform(other.coeffRef(0,0), other.coeffRef(1,0), other.coeffRef(2,0)
other.coeffRef(0,1), other.coeffRef(1,1), other.coeffRef(2,1)
other.coeffRef(0,2), other.coeffRef(1,2), other.coeffRef(2,2);
@@ -375,7 +375,7 @@ template<typename OtherDerived>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::scale(const MatrixBase<OtherDerived> &other)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim));
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
linear() = (linear() * other.asDiagonal()).lazy();
return *this;
}
@@ -400,7 +400,7 @@ template<typename OtherDerived>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::prescale(const MatrixBase<OtherDerived> &other)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim));
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
m_matrix.template block<Dim,HDim>(0,0) = (other.asDiagonal() * m_matrix.template block<Dim,HDim>(0,0)).lazy();
return *this;
}
@@ -425,7 +425,7 @@ template<typename OtherDerived>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::translate(const MatrixBase<OtherDerived> &other)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim));
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
translation() += linear() * other;
return *this;
}
@@ -439,7 +439,7 @@ template<typename OtherDerived>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::pretranslate(const MatrixBase<OtherDerived> &other)
{
- EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim));
+ EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
translation() += other;
return *this;
}
@@ -496,7 +496,7 @@ template<typename Scalar, int Dim>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::shear(Scalar sx, Scalar sy)
{
- EIGEN_STATIC_ASSERT(int(Dim)==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(int(Dim)==2, you_made_a_programming_mistake)
VectorType tmp = linear().col(0)*sy + linear().col(1);
linear() << linear().col(0) + linear().col(1)*sx, tmp;
return *this;
@@ -511,7 +511,7 @@ template<typename Scalar, int Dim>
Transform<Scalar,Dim>&
Transform<Scalar,Dim>::preshear(Scalar sx, Scalar sy)
{
- EIGEN_STATIC_ASSERT(int(Dim)==2, you_did_a_programming_error);
+ EIGEN_STATIC_ASSERT(int(Dim)==2, you_made_a_programming_mistake)
m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0);
return *this;
}
diff --git a/Eigen/src/LU/Inverse.h b/Eigen/src/LU/Inverse.h
index e08ece2c4..cd59290ab 100644
--- a/Eigen/src/LU/Inverse.h
+++ b/Eigen/src/LU/Inverse.h
@@ -220,7 +220,7 @@ inline void MatrixBase<Derived>::computeInverse(EvalType *result) const
{
typedef typename ei_eval<Derived>::type MatrixType;
ei_assert(rows() == cols());
- EIGEN_STATIC_ASSERT(NumTraits<Scalar>::HasFloatingPoint,scalar_type_must_be_floating_point);
+ EIGEN_STATIC_ASSERT(NumTraits<Scalar>::HasFloatingPoint,numeric_type_must_be_floating_point)
ei_compute_inverse<MatrixType>::run(eval(), result);
}
diff --git a/Eigen/src/Sparse/SparseBlock.h b/Eigen/src/Sparse/SparseBlock.h
index dd952714d..3b1abf13f 100644
--- a/Eigen/src/Sparse/SparseBlock.h
+++ b/Eigen/src/Sparse/SparseBlock.h
@@ -59,7 +59,7 @@ public:
: m_matrix(matrix), m_startRow(startRow), m_startCol(startCol),
m_blockRows(matrix.rows()), m_blockCols(matrix.cols())
{
- EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic,this_method_is_only_for_fixed_size);
+ EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic,this_method_is_only_for_fixed_size)
ei_assert(startRow >= 0 && BlockRows >= 1 && startRow + BlockRows <= matrix.rows()
&& startCol >= 0 && BlockCols >= 1 && startCol + BlockCols <= matrix.cols());
}
diff --git a/test/regression.cpp b/test/regression.cpp
index d9197bbfc..5893c9564 100644
--- a/test/regression.cpp
+++ b/test/regression.cpp
@@ -53,8 +53,8 @@ void makeNoisyCohyperplanarPoints(int numPoints,
// project cur_point onto the hyperplane
Scalar x = - (hyperplane->coeffs().start(size).cwise()*cur_point).sum();
cur_point *= hyperplane->coeffs().coeff(size) / x;
- } while( ei_abs(cur_point.norm()) < 0.5
- || ei_abs(cur_point.norm()) > 2.0 );
+ } while( cur_point.norm() < 0.5
+ || cur_point.norm() > 2.0 );
}
// add some noise to these points
@@ -96,7 +96,7 @@ void test_regression()
Vector4d points4d [1000];
Vector4d *points4d_ptrs [1000];
for(int i = 0; i < 1000; i++) points4d_ptrs[i] = &(points4d[i]);
- Hyperplane<float,4> coeffs5d;
+ Hyperplane<double,4> coeffs5d;
makeNoisyCohyperplanarPoints(1000, points4d_ptrs, &coeffs5d, 0.01);
CALL_SUBTEST(check_fitHyperplane(10, points4d_ptrs, coeffs5d, 0.05));
CALL_SUBTEST(check_fitHyperplane(100, points4d_ptrs, coeffs5d, 0.01));