aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Array/Replicate.h2
-rw-r--r--Eigen/src/Array/Reverse.h4
-rw-r--r--Eigen/src/Array/Select.h3
-rw-r--r--Eigen/src/Array/VectorwiseOp.h3
-rw-r--r--Eigen/src/Cholesky/LDLT.h2
-rw-r--r--Eigen/src/Core/Block.h4
-rw-r--r--Eigen/src/Core/Coeffs.h2
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h1
-rw-r--r--Eigen/src/Core/CwiseNullaryOp.h3
-rw-r--r--Eigen/src/Core/CwiseUnaryOp.h3
-rw-r--r--Eigen/src/Core/CwiseUnaryView.h3
-rw-r--r--Eigen/src/Core/Diagonal.h3
-rw-r--r--Eigen/src/Core/DiagonalProduct.h3
-rw-r--r--Eigen/src/Core/Minor.h3
-rw-r--r--Eigen/src/Core/ProductBase.h5
-rw-r--r--Eigen/src/Core/Redux.h2
-rw-r--r--Eigen/src/Core/SelfAdjointView.h5
-rw-r--r--Eigen/src/Core/SelfCwiseBinaryOp.h2
-rw-r--r--Eigen/src/Core/StableNorm.h4
-rw-r--r--Eigen/src/Core/Transpose.h3
-rw-r--r--Eigen/src/Core/TriangularMatrix.h3
-rw-r--r--Eigen/src/Core/products/GeneralUnrolled.h4
-rw-r--r--Eigen/src/Core/util/Constants.h7
-rw-r--r--Eigen/src/Core/util/XprHelper.h53
-rw-r--r--Eigen/src/Geometry/Homogeneous.h3
-rw-r--r--Eigen/src/LU/FullPivLU.h2
-rw-r--r--Eigen/src/Sparse/DynamicSparseMatrix.h9
-rw-r--r--Eigen/src/Sparse/SparseMatrix.h9
-rw-r--r--Eigen/src/Sparse/SparseVector.h9
29 files changed, 84 insertions, 75 deletions
diff --git a/Eigen/src/Array/Replicate.h b/Eigen/src/Array/Replicate.h
index 0b1b58f1b..7e7595e15 100644
--- a/Eigen/src/Array/Replicate.h
+++ b/Eigen/src/Array/Replicate.h
@@ -55,7 +55,7 @@ struct ei_traits<Replicate<MatrixType,RowFactor,ColFactor> >
: ColFactor * MatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
- Flags = _MatrixTypeNested::Flags & HereditaryBits,
+ Flags = (_MatrixTypeNested::Flags & HereditaryBits) | EIGEN_PROPAGATE_NESTING_BIT(MatrixType::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/Array/Reverse.h b/Eigen/src/Array/Reverse.h
index c908e8d9c..24fcf644b 100644
--- a/Eigen/src/Array/Reverse.h
+++ b/Eigen/src/Array/Reverse.h
@@ -59,7 +59,9 @@ struct ei_traits<Reverse<MatrixType, Direction> >
LinearAccess = ( (Direction==BothDirections) && (int(_MatrixTypeNested::Flags)&PacketAccessBit) )
? LinearAccessBit : 0,
- Flags = (int(_MatrixTypeNested::Flags) & (HereditaryBits | PacketAccessBit | LinearAccess)),
+ Flags = (int(_MatrixTypeNested::Flags)
+ & (HereditaryBits | PacketAccessBit | LinearAccess))
+ | EIGEN_PROPAGATE_NESTING_BIT(MatrixType::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
diff --git a/Eigen/src/Array/Select.h b/Eigen/src/Array/Select.h
index 38c1a716f..bcb233c80 100644
--- a/Eigen/src/Array/Select.h
+++ b/Eigen/src/Array/Select.h
@@ -55,7 +55,8 @@ struct ei_traits<Select<ConditionMatrixType, ThenMatrixType, ElseMatrixType> >
ColsAtCompileTime = ConditionMatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = ConditionMatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = ConditionMatrixType::MaxColsAtCompileTime,
- Flags = (unsigned int)ThenMatrixType::Flags & ElseMatrixType::Flags & HereditaryBits,
+ Flags = ((unsigned int)ThenMatrixType::Flags & ElseMatrixType::Flags & HereditaryBits)
+ | EIGEN_PROPAGATE_NESTING_BIT(ConditionMatrixType::Flags|ThenMatrixType::Flags|ElseMatrixType::Flags),
CoeffReadCost = ei_traits<typename ei_cleantype<ConditionMatrixNested>::type>::CoeffReadCost
+ EIGEN_ENUM_MAX(ei_traits<typename ei_cleantype<ThenMatrixNested>::type>::CoeffReadCost,
ei_traits<typename ei_cleantype<ElseMatrixNested>::type>::CoeffReadCost)
diff --git a/Eigen/src/Array/VectorwiseOp.h b/Eigen/src/Array/VectorwiseOp.h
index 697a07d32..bf5750b13 100644
--- a/Eigen/src/Array/VectorwiseOp.h
+++ b/Eigen/src/Array/VectorwiseOp.h
@@ -60,7 +60,8 @@ struct ei_traits<PartialReduxExpr<MatrixType, MemberOp, Direction> >
ColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = Direction==Vertical ? 1 : MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = Direction==Horizontal ? 1 : MatrixType::MaxColsAtCompileTime,
- Flags = (unsigned int)_MatrixTypeNested::Flags & HereditaryBits,
+ Flags = ((unsigned int)_MatrixTypeNested::Flags & HereditaryBits)
+ | EIGEN_PROPAGATE_NESTING_BIT(_MatrixTypeNested::Flags),
TraversalSize = Direction==Vertical ? RowsAtCompileTime : ColsAtCompileTime
};
#if EIGEN_GNUC_AT_LEAST(3,4)
diff --git a/Eigen/src/Cholesky/LDLT.h b/Eigen/src/Cholesky/LDLT.h
index 8b13d8e2e..54bd2d58a 100644
--- a/Eigen/src/Cholesky/LDLT.h
+++ b/Eigen/src/Cholesky/LDLT.h
@@ -206,7 +206,7 @@ LDLT<MatrixType>& LDLT<MatrixType>::compute(const MatrixType& a)
// in "Analysis of the Cholesky Decomposition of a Semi-definite Matrix" by
// Nicholas J. Higham. Also see "Accuracy and Stability of Numerical
// Algorithms" page 217, also by Higham.
- cutoff = ei_abs(epsilon<Scalar>() * size * biggest_in_corner);
+ cutoff = ei_abs(epsilon<Scalar>() * RealScalar(size) * biggest_in_corner);
m_sign = ei_real(m_matrix.diagonal().coeff(index_of_biggest_in_corner)) > 0 ? 1 : -1;
}
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index 3b4234c22..244f33ca1 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -76,7 +76,9 @@ struct ei_traits<Block<MatrixType, BlockRows, BlockCols, _DirectAccessStatus> >
MaskPacketAccessBit = (InnerMaxSize == Dynamic || (InnerSize >= ei_packet_traits<Scalar>::size))
? PacketAccessBit : 0,
FlagsLinearAccessBit = (RowsAtCompileTime == 1 || ColsAtCompileTime == 1) ? LinearAccessBit : 0,
- Flags = (ei_traits<MatrixType>::Flags & (HereditaryBits | MaskPacketAccessBit | DirectAccessBit)) | FlagsLinearAccessBit
+ Flags = (ei_traits<MatrixType>::Flags & (HereditaryBits | MaskPacketAccessBit | DirectAccessBit))
+ | FlagsLinearAccessBit
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags)
};
};
diff --git a/Eigen/src/Core/Coeffs.h b/Eigen/src/Core/Coeffs.h
index ebfd0c80e..77f0ea544 100644
--- a/Eigen/src/Core/Coeffs.h
+++ b/Eigen/src/Core/Coeffs.h
@@ -343,7 +343,7 @@ template<typename OtherDerived>
EIGEN_STRONG_INLINE void DenseBase<Derived>::copyCoeff(int index, const DenseBase<OtherDerived>& other)
{
ei_internal_assert(index >= 0 && index < size());
- derived().coeffRef(index) = other.derived().coeff(index);
+ derived().coeffRef(index) = Scalar(other.derived().coeff(index));
}
/** \internal Copies the packet at position (row,col) of other into *this.
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index 9ed005dce..a39f206cf 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -72,6 +72,7 @@ struct ei_traits<CwiseBinaryOp<BinaryOp, Lhs, Rhs> > : ei_traits<Lhs>
( AlignedBit
| (StorageOrdersAgree ? LinearAccessBit : 0)
| (ei_functor_traits<BinaryOp>::PacketAccess && StorageOrdersAgree ? PacketAccessBit : 0)
+ | EIGEN_PROPAGATE_NESTING_BIT(Lhs::Flags|Rhs::Flags)
)
)
),
diff --git a/Eigen/src/Core/CwiseNullaryOp.h b/Eigen/src/Core/CwiseNullaryOp.h
index 5800335d7..e4880e5d0 100644
--- a/Eigen/src/Core/CwiseNullaryOp.h
+++ b/Eigen/src/Core/CwiseNullaryOp.h
@@ -48,7 +48,8 @@ struct ei_traits<CwiseNullaryOp<NullaryOp, MatrixType> > : ei_traits<MatrixType>
& ( HereditaryBits
| (ei_functor_has_linear_access<NullaryOp>::ret ? LinearAccessBit : 0)
| (ei_functor_traits<NullaryOp>::PacketAccess ? PacketAccessBit : 0)))
- | (ei_functor_traits<NullaryOp>::IsRepeatable ? 0 : EvalBeforeNestingBit),
+ | (ei_functor_traits<NullaryOp>::IsRepeatable ? 0 : EvalBeforeNestingBit)
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = ei_functor_traits<NullaryOp>::Cost
};
};
diff --git a/Eigen/src/Core/CwiseUnaryOp.h b/Eigen/src/Core/CwiseUnaryOp.h
index b51bd51af..b519e3b1d 100644
--- a/Eigen/src/Core/CwiseUnaryOp.h
+++ b/Eigen/src/Core/CwiseUnaryOp.h
@@ -51,7 +51,8 @@ struct ei_traits<CwiseUnaryOp<UnaryOp, MatrixType> >
enum {
Flags = (_MatrixTypeNested::Flags & (
HereditaryBits | LinearAccessBit | AlignedBit
- | (ei_functor_traits<UnaryOp>::PacketAccess ? PacketAccessBit : 0))),
+ | (ei_functor_traits<UnaryOp>::PacketAccess ? PacketAccessBit : 0)))
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost + ei_functor_traits<UnaryOp>::Cost
};
};
diff --git a/Eigen/src/Core/CwiseUnaryView.h b/Eigen/src/Core/CwiseUnaryView.h
index 2198ed226..deef14597 100644
--- a/Eigen/src/Core/CwiseUnaryView.h
+++ b/Eigen/src/Core/CwiseUnaryView.h
@@ -47,7 +47,8 @@ struct ei_traits<CwiseUnaryView<ViewOp, MatrixType> >
typedef typename MatrixType::Nested MatrixTypeNested;
typedef typename ei_cleantype<MatrixTypeNested>::type _MatrixTypeNested;
enum {
- Flags = (ei_traits<_MatrixTypeNested>::Flags & (HereditaryBits | LinearAccessBit | AlignedBit)),
+ Flags = (ei_traits<_MatrixTypeNested>::Flags & (HereditaryBits | LinearAccessBit | AlignedBit))
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags), // if I am not wrong, I need to test this on MatrixType and not on the nested type
CoeffReadCost = ei_traits<_MatrixTypeNested>::CoeffReadCost + ei_functor_traits<ViewOp>::Cost
};
};
diff --git a/Eigen/src/Core/Diagonal.h b/Eigen/src/Core/Diagonal.h
index 3720952cd..4f8bc01c3 100644
--- a/Eigen/src/Core/Diagonal.h
+++ b/Eigen/src/Core/Diagonal.h
@@ -58,7 +58,8 @@ struct ei_traits<Diagonal<MatrixType,Index> >
: Index == Dynamic ? EIGEN_ENUM_MIN(MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime)
: (EIGEN_ENUM_MIN(MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime) - AbsIndex),
MaxColsAtCompileTime = 1,
- Flags = (unsigned int)_MatrixTypeNested::Flags & (HereditaryBits | LinearAccessBit),
+ Flags = ((unsigned int)_MatrixTypeNested::Flags & (HereditaryBits | LinearAccessBit))
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/Core/DiagonalProduct.h b/Eigen/src/Core/DiagonalProduct.h
index 868b4419a..7871e7704 100644
--- a/Eigen/src/Core/DiagonalProduct.h
+++ b/Eigen/src/Core/DiagonalProduct.h
@@ -37,7 +37,8 @@ struct ei_traits<DiagonalProduct<MatrixType, DiagonalType, ProductOrder> >
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
Flags = (HereditaryBits & (unsigned int)(MatrixType::Flags))
- | (PacketAccessBit & (unsigned int)(MatrixType::Flags) & (unsigned int)(DiagonalType::DiagonalVectorType::Flags)),
+ | (PacketAccessBit & (unsigned int)(MatrixType::Flags) & (unsigned int)(DiagonalType::DiagonalVectorType::Flags))
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = NumTraits<Scalar>::MulCost + MatrixType::CoeffReadCost + DiagonalType::DiagonalVectorType::CoeffReadCost
};
};
diff --git a/Eigen/src/Core/Minor.h b/Eigen/src/Core/Minor.h
index e7e164a16..aa7739370 100644
--- a/Eigen/src/Core/Minor.h
+++ b/Eigen/src/Core/Minor.h
@@ -53,7 +53,8 @@ struct ei_traits<Minor<MatrixType> >
int(MatrixType::MaxRowsAtCompileTime) - 1 : Dynamic,
MaxColsAtCompileTime = (MatrixType::MaxColsAtCompileTime != Dynamic) ?
int(MatrixType::MaxColsAtCompileTime) - 1 : Dynamic,
- Flags = _MatrixTypeNested::Flags & HereditaryBits,
+ Flags = (_MatrixTypeNested::Flags & HereditaryBits)
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost // minor is used typically on tiny matrices,
// where loops are unrolled and the 'if' evaluates at compile time
};
diff --git a/Eigen/src/Core/ProductBase.h b/Eigen/src/Core/ProductBase.h
index 5c51ea27c..94f01f55e 100644
--- a/Eigen/src/Core/ProductBase.h
+++ b/Eigen/src/Core/ProductBase.h
@@ -42,7 +42,10 @@ struct ei_traits<ProductBase<Derived,_Lhs,_Rhs> > //: ei_traits<typename ei_clea
ColsAtCompileTime = ei_traits<Rhs>::ColsAtCompileTime,
MaxRowsAtCompileTime = ei_traits<Lhs>::MaxRowsAtCompileTime,
MaxColsAtCompileTime = ei_traits<Rhs>::MaxColsAtCompileTime,
- Flags = EvalBeforeNestingBit | EvalBeforeAssigningBit,
+ Flags = EvalBeforeNestingBit
+ | EvalBeforeAssigningBit
+ | NestParentByRefBit
+ | EIGEN_PROPAGATE_NESTING_BIT(Lhs::Flags|Rhs::Flags),
CoeffReadCost = 0 // FIXME why is it needed ?
};
};
diff --git a/Eigen/src/Core/Redux.h b/Eigen/src/Core/Redux.h
index 42287dab9..d5b0c60c2 100644
--- a/Eigen/src/Core/Redux.h
+++ b/Eigen/src/Core/Redux.h
@@ -355,7 +355,7 @@ template<typename Derived>
EIGEN_STRONG_INLINE typename ei_traits<Derived>::Scalar
DenseBase<Derived>::mean() const
{
- return this->redux(Eigen::ei_scalar_sum_op<Scalar>()) / Scalar(this->size());
+ return Scalar(this->redux(Eigen::ei_scalar_sum_op<Scalar>())) / Scalar(this->size());
}
/** \returns the product of all coefficients of *this
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index 6d01ee495..9f1893c4a 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -47,8 +47,9 @@ struct ei_traits<SelfAdjointView<MatrixType, TriangularPart> > : ei_traits<Matri
typedef MatrixType ExpressionType;
enum {
Mode = TriangularPart | SelfAdjoint,
- Flags = _MatrixTypeNested::Flags & (HereditaryBits)
- & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)), // FIXME these flags should be preserved
+ Flags = (_MatrixTypeNested::Flags & (HereditaryBits)
+ & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) // FIXME these flags should be preserved
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/Core/SelfCwiseBinaryOp.h b/Eigen/src/Core/SelfCwiseBinaryOp.h
index 7ae2e82a4..c403a4f55 100644
--- a/Eigen/src/Core/SelfCwiseBinaryOp.h
+++ b/Eigen/src/Core/SelfCwiseBinaryOp.h
@@ -90,7 +90,7 @@ template<typename BinaryOp, typename MatrixType> class SelfCwiseBinaryOp
OtherDerived& _other = other.const_cast_derived();
ei_internal_assert(index >= 0 && index < m_matrix.size());
Scalar& tmp = m_matrix.coeffRef(index);
- tmp = m_functor(tmp, _other.coeff(index));
+ tmp = m_functor(tmp, Scalar(_other.coeff(index)));
}
template<typename OtherDerived, int StoreMode, int LoadMode>
diff --git a/Eigen/src/Core/StableNorm.h b/Eigen/src/Core/StableNorm.h
index b4d6aa353..70371d8ba 100644
--- a/Eigen/src/Core/StableNorm.h
+++ b/Eigen/src/Core/StableNorm.h
@@ -118,10 +118,10 @@ MatrixBase<Derived>::blueNorm() const
iexp = - ((iemax+it)/2);
s2m = RealScalar(std::pow(RealScalar(ibeta),RealScalar(iexp))); // scaling factor for upper range
- overfl = rbig*s2m; // overfow boundary for abig
+ overfl = rbig*s2m; // overflow boundary for abig
eps = RealScalar(std::pow(double(ibeta), 1-it));
relerr = ei_sqrt(eps); // tolerance for neglecting asml
- abig = 1.0/eps - 1.0;
+ abig = RealScalar(1.0/eps - 1.0);
if (RealScalar(nbig)>abig) nmax = int(abig); // largest safe n
else nmax = nbig;
}
diff --git a/Eigen/src/Core/Transpose.h b/Eigen/src/Core/Transpose.h
index 143b39033..18e4a1739 100644
--- a/Eigen/src/Core/Transpose.h
+++ b/Eigen/src/Core/Transpose.h
@@ -50,7 +50,8 @@ struct ei_traits<Transpose<MatrixType> > : ei_traits<MatrixType>
ColsAtCompileTime = MatrixType::RowsAtCompileTime,
MaxRowsAtCompileTime = MatrixType::MaxColsAtCompileTime,
MaxColsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
- Flags = (int(_MatrixTypeNested::Flags) ^ RowMajorBit),
+ Flags = (int(_MatrixTypeNested::Flags & ~NestByRefBit) ^ RowMajorBit)
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index 8bea0aa68..e6690b6fa 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -130,7 +130,8 @@ struct ei_traits<TriangularView<MatrixType, _Mode> > : ei_traits<MatrixType>
typedef MatrixType ExpressionType;
enum {
Mode = _Mode,
- Flags = (_MatrixTypeNested::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) | Mode,
+ Flags = (_MatrixTypeNested::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) | Mode
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/Core/products/GeneralUnrolled.h b/Eigen/src/Core/products/GeneralUnrolled.h
index 32aa3afe6..a8fec1fac 100644
--- a/Eigen/src/Core/products/GeneralUnrolled.h
+++ b/Eigen/src/Core/products/GeneralUnrolled.h
@@ -81,8 +81,10 @@ struct ei_traits<GeneralProduct<LhsNested,RhsNested,UnrolledProduct> >
Flags = ((unsigned int)(LhsFlags | RhsFlags) & HereditaryBits & RemovedBits)
| EvalBeforeAssigningBit
| EvalBeforeNestingBit
+ | NestParentByRefBit
| (CanVectorizeLhs || CanVectorizeRhs ? PacketAccessBit : 0)
- | (LhsFlags & RhsFlags & AlignedBit),
+ | (LhsFlags & RhsFlags & AlignedBit)
+ | EIGEN_PROPAGATE_NESTING_BIT(LhsFlags|RhsFlags),
CoeffReadCost = InnerSize == Dynamic ? Dynamic
: InnerSize * (NumTraits<Scalar>::MulCost + LhsCoeffReadCost + RhsCoeffReadCost)
diff --git a/Eigen/src/Core/util/Constants.h b/Eigen/src/Core/util/Constants.h
index c747d970b..1b34b401d 100644
--- a/Eigen/src/Core/util/Constants.h
+++ b/Eigen/src/Core/util/Constants.h
@@ -35,7 +35,7 @@
* - 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.
+ * it would be a meta-selector. Doing this everywhere would reduce code readability and lengthen compilation times.
* Also, disabling compiler warnings for integer overflow, sounds like a bad idea.
* - It should be a prime number, because for example the old value 10000 led to bugs with 100x100 matrices.
*
@@ -76,7 +76,7 @@ const unsigned int EvalBeforeNestingBit = 0x2;
/** \ingroup flags
*
- * means the expression should be evaluated before any assignement */
+ * means the expression should be evaluated before any assignment */
const unsigned int EvalBeforeAssigningBit = 0x4;
/** \ingroup flags
@@ -97,6 +97,9 @@ const unsigned int EvalBeforeAssigningBit = 0x4;
*/
const unsigned int PacketAccessBit = 0x8;
+const unsigned int NestParentByRefBit = 0x80;
+const unsigned int NestByRefBit = 0x100;
+
#ifdef EIGEN_VECTORIZE
/** \ingroup flags
*
diff --git a/Eigen/src/Core/util/XprHelper.h b/Eigen/src/Core/util/XprHelper.h
index 136cc876b..77b3968b1 100644
--- a/Eigen/src/Core/util/XprHelper.h
+++ b/Eigen/src/Core/util/XprHelper.h
@@ -97,7 +97,7 @@ class ei_compute_matrix_flags
};
public:
- enum { ret = LinearAccessBit | DirectAccessBit | packet_access_bit | row_major_bit | aligned_bit };
+ enum { ret = LinearAccessBit | DirectAccessBit | NestByRefBit | packet_access_bit | row_major_bit | aligned_bit };
};
template<int _Rows, int _Cols> struct ei_size_at_compile_time
@@ -201,6 +201,28 @@ template<typename T> struct ei_plain_matrix_type_row_major
// we should be able to get rid of this one too
template<typename T> struct ei_must_nest_by_value { enum { ret = false }; };
+template<class T>
+struct ei_is_reference
+{
+#ifndef NDEBUG
+ static void check() { std::cout << typeid(T).name() << std::endl; }
+#else
+ static void check() {}
+#endif
+ enum { ret = false };
+};
+
+template<class T>
+struct ei_is_reference<T&>
+{
+#ifndef NDEBUG
+ static void check() { std::cout << typeid(T).name() << "&" << std::endl; }
+#else
+ static void check() {}
+#endif
+ enum { ret = true };
+};
+
/**
* The reference selector for template expressions. The idea is that we don't
* need to use references for expressions since they are light weight proxy
@@ -209,31 +231,14 @@ template<typename T> struct ei_must_nest_by_value { enum { ret = false }; };
template <typename T>
struct ei_ref_selector
{
- typedef T type;
-};
-
-/**
-* Matrices on the other hand side should only be copied, when it is sure
-* we gain by copying (see arithmetic cost check and eval before nesting flag).
-* Note: This is an optimization measure that comprises potential (though little)
-* to create erroneous code. Any user, utilizing ei_nested outside of
-* Eigen needs to take care that no references to temporaries are
-* stored or that this potential danger is at least communicated
-* to the user.
-**/
-template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
-struct ei_ref_selector< Matrix<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> >
-{
- typedef Matrix<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> MatrixType;
- typedef MatrixType const& type;
+ typedef typename ei_meta_if<
+ bool(ei_traits<T>::Flags & NestByRefBit),
+ T const&,
+ T
+ >::ret type;
};
-template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
-struct ei_ref_selector< Array<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> >
-{
- typedef Array<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> ArrayType;
- typedef ArrayType const& type;
-};
+#define EIGEN_PROPAGATE_NESTING_BIT(ReferenceFlags) ((ReferenceFlags) & NestParentByRefBit)<<1
/** \internal Determines how a given expression should be nested into another one.
* For example, when you do a * (b+c), Eigen will determine how the expression b+c should be
diff --git a/Eigen/src/Geometry/Homogeneous.h b/Eigen/src/Geometry/Homogeneous.h
index 76ca66c57..8c95a4d4d 100644
--- a/Eigen/src/Geometry/Homogeneous.h
+++ b/Eigen/src/Geometry/Homogeneous.h
@@ -55,7 +55,8 @@ struct ei_traits<Homogeneous<MatrixType,Direction> >
ColsAtCompileTime = Direction==Horizontal ? ColsPlusOne : MatrixType::ColsAtCompileTime,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
- Flags = _MatrixTypeNested::Flags & HereditaryBits,
+ Flags = _MatrixTypeNested::Flags & HereditaryBits
+ | EIGEN_PROPAGATE_NESTING_BIT(ei_traits<MatrixType>::Flags),
CoeffReadCost = _MatrixTypeNested::CoeffReadCost
};
};
diff --git a/Eigen/src/LU/FullPivLU.h b/Eigen/src/LU/FullPivLU.h
index fe14a9080..657cbbf1c 100644
--- a/Eigen/src/LU/FullPivLU.h
+++ b/Eigen/src/LU/FullPivLU.h
@@ -476,7 +476,7 @@ typename ei_traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() cons
{
ei_assert(m_isInitialized && "LU is not initialized.");
ei_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
- return Scalar(m_det_pq) * m_lu.diagonal().prod();
+ return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
}
/********* Implementation of kernel() **************************************************/
diff --git a/Eigen/src/Sparse/DynamicSparseMatrix.h b/Eigen/src/Sparse/DynamicSparseMatrix.h
index 4373e1cda..96385df50 100644
--- a/Eigen/src/Sparse/DynamicSparseMatrix.h
+++ b/Eigen/src/Sparse/DynamicSparseMatrix.h
@@ -52,19 +52,12 @@ struct ei_traits<DynamicSparseMatrix<_Scalar, _Flags> >
ColsAtCompileTime = Dynamic,
MaxRowsAtCompileTime = Dynamic,
MaxColsAtCompileTime = Dynamic,
- Flags = _Flags,
+ Flags = _Flags | NestByRefBit,
CoeffReadCost = NumTraits<Scalar>::ReadCost,
SupportedAccessPatterns = OuterRandomAccessPattern
};
};
-template<typename _Scalar, int _Options>
-struct ei_ref_selector< DynamicSparseMatrix<_Scalar, _Options> >
-{
- typedef DynamicSparseMatrix<_Scalar, _Options> MatrixType;
- typedef MatrixType const& type;
-};
-
template<typename _Scalar, int _Flags>
class DynamicSparseMatrix
: public SparseMatrixBase<DynamicSparseMatrix<_Scalar, _Flags> >
diff --git a/Eigen/src/Sparse/SparseMatrix.h b/Eigen/src/Sparse/SparseMatrix.h
index d884e7df8..bba0404c3 100644
--- a/Eigen/src/Sparse/SparseMatrix.h
+++ b/Eigen/src/Sparse/SparseMatrix.h
@@ -51,20 +51,13 @@ struct ei_traits<SparseMatrix<_Scalar, _Options> >
ColsAtCompileTime = Dynamic,
MaxRowsAtCompileTime = Dynamic,
MaxColsAtCompileTime = Dynamic,
- Flags = _Options,
+ Flags = _Options | NestByRefBit,
CoeffReadCost = NumTraits<Scalar>::ReadCost,
SupportedAccessPatterns = InnerRandomAccessPattern
};
};
template<typename _Scalar, int _Options>
-struct ei_ref_selector<SparseMatrix<_Scalar, _Options> >
-{
- typedef SparseMatrix<_Scalar, _Options> MatrixType;
- typedef MatrixType const& type;
-};
-
-template<typename _Scalar, int _Options>
class SparseMatrix
: public SparseMatrixBase<SparseMatrix<_Scalar, _Options> >
{
diff --git a/Eigen/src/Sparse/SparseVector.h b/Eigen/src/Sparse/SparseVector.h
index 5ac3f6be7..13b4f75a6 100644
--- a/Eigen/src/Sparse/SparseVector.h
+++ b/Eigen/src/Sparse/SparseVector.h
@@ -46,20 +46,13 @@ struct ei_traits<SparseVector<_Scalar, _Options> >
ColsAtCompileTime = IsColVector ? 1 : Dynamic,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
- Flags = _Options,
+ Flags = _Options | NestByRefBit,
CoeffReadCost = NumTraits<Scalar>::ReadCost,
SupportedAccessPatterns = InnerRandomAccessPattern
};
};
template<typename _Scalar, int _Options>
-struct ei_ref_selector< SparseVector<_Scalar, _Options> >
-{
- typedef SparseVector<_Scalar, _Options> MatrixType;
- typedef MatrixType const& type;
-};
-
-template<typename _Scalar, int _Options>
class SparseVector
: public SparseMatrixBase<SparseVector<_Scalar, _Options> >
{