aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-06-09 09:11:12 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-06-09 09:11:12 +0200
commit9a2447b0c95bb62bd3e5ae68bf14d2c5e66a18cf (patch)
tree40712813d830e39c1e8e5b32be7c46abf4c90695
parentcd8b996f99de67035a0504cbaf0a627fb68f0f1d (diff)
Fix shadow warnings triggered by clang
-rw-r--r--Eigen/src/Core/arch/AVX/PacketMath.h20
-rw-r--r--Eigen/src/Core/products/TriangularSolverMatrix.h4
-rw-r--r--Eigen/src/Eigenvalues/RealQZ.h6
-rw-r--r--Eigen/src/SVD/BDCSVD.h2
-rw-r--r--Eigen/src/SVD/JacobiSVD.h2
-rw-r--r--Eigen/src/SparseCore/SparseMatrix.h12
-rw-r--r--Eigen/src/SparseCore/SparseView.h6
-rw-r--r--test/inverse.cpp10
-rw-r--r--test/sparse_solver.h1
-rw-r--r--test/sparse_vector.cpp10
10 files changed, 38 insertions, 35 deletions
diff --git a/Eigen/src/Core/arch/AVX/PacketMath.h b/Eigen/src/Core/arch/AVX/PacketMath.h
index 1e751dc32..a2306fd1a 100644
--- a/Eigen/src/Core/arch/AVX/PacketMath.h
+++ b/Eigen/src/Core/arch/AVX/PacketMath.h
@@ -434,26 +434,30 @@ struct palign_impl<Offset,Packet8f>
if (Offset==1)
{
first = _mm256_blend_ps(first, second, 1);
- Packet8f tmp = _mm256_permute_ps (first, _MM_SHUFFLE(0,3,2,1));
- first = _mm256_blend_ps(tmp, _mm256_permute2f128_ps (tmp, tmp, 1), 0x88);
+ Packet8f tmp1 = _mm256_permute_ps (first, _MM_SHUFFLE(0,3,2,1));
+ Packet8f tmp2 = _mm256_permute2f128_ps (tmp1, tmp1, 1);
+ first = _mm256_blend_ps(tmp1, tmp2, 0x88);
}
else if (Offset==2)
{
first = _mm256_blend_ps(first, second, 3);
- Packet8f tmp = _mm256_permute_ps (first, _MM_SHUFFLE(1,0,3,2));
- first = _mm256_blend_ps(tmp, _mm256_permute2f128_ps (tmp, tmp, 1), 0xcc);
+ Packet8f tmp1 = _mm256_permute_ps (first, _MM_SHUFFLE(1,0,3,2));
+ Packet8f tmp2 = _mm256_permute2f128_ps (tmp1, tmp1, 1);
+ first = _mm256_blend_ps(tmp1, tmp2, 0xcc);
}
else if (Offset==3)
{
first = _mm256_blend_ps(first, second, 7);
- Packet8f tmp = _mm256_permute_ps (first, _MM_SHUFFLE(2,1,0,3));
- first = _mm256_blend_ps(tmp, _mm256_permute2f128_ps (tmp, tmp, 1), 0xee);
+ Packet8f tmp1 = _mm256_permute_ps (first, _MM_SHUFFLE(2,1,0,3));
+ Packet8f tmp2 = _mm256_permute2f128_ps (tmp1, tmp1, 1);
+ first = _mm256_blend_ps(tmp1, tmp2, 0xee);
}
else if (Offset==4)
{
first = _mm256_blend_ps(first, second, 15);
- Packet8f tmp = _mm256_permute_ps (first, _MM_SHUFFLE(3,2,1,0));
- first = _mm256_permute_ps(_mm256_permute2f128_ps (tmp, tmp, 1), _MM_SHUFFLE(3,2,1,0));
+ Packet8f tmp1 = _mm256_permute_ps (first, _MM_SHUFFLE(3,2,1,0));
+ Packet8f tmp2 = _mm256_permute2f128_ps (tmp1, tmp1, 1);
+ first = _mm256_permute_ps(tmp2, _MM_SHUFFLE(3,2,1,0));
}
else if (Offset==5)
{
diff --git a/Eigen/src/Core/products/TriangularSolverMatrix.h b/Eigen/src/Core/products/TriangularSolverMatrix.h
index f5de67c59..a9a198d64 100644
--- a/Eigen/src/Core/products/TriangularSolverMatrix.h
+++ b/Eigen/src/Core/products/TriangularSolverMatrix.h
@@ -117,8 +117,9 @@ EIGEN_DONT_INLINE void triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conju
{
// TODO write a small kernel handling this (can be shared with trsv)
Index i = IsLower ? k2+k1+k : k2-k1-k-1;
- Index s = IsLower ? k2+k1 : i+1;
Index rs = actualPanelWidth - k - 1; // remaining size
+ Index s = TriStorageOrder==RowMajor ? (IsLower ? k2+k1 : i+1)
+ : IsLower ? i+1 : i-rs;
Scalar a = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(tri(i,i));
for (Index j=j2; j<j2+actual_cols; ++j)
@@ -135,7 +136,6 @@ EIGEN_DONT_INLINE void triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conju
}
else
{
- Index s = IsLower ? i+1 : i-rs;
Scalar b = (other(i,j) *= a);
Scalar* r = &other(s,j);
const Scalar* l = &tri(s,i);
diff --git a/Eigen/src/Eigenvalues/RealQZ.h b/Eigen/src/Eigenvalues/RealQZ.h
index 677c7c0bb..02ebb7d17 100644
--- a/Eigen/src/Eigenvalues/RealQZ.h
+++ b/Eigen/src/Eigenvalues/RealQZ.h
@@ -315,8 +315,8 @@ namespace Eigen {
const Index dim=m_S.cols();
if (abs(m_S.coeff(i+1,i))==Scalar(0))
return;
- Index z = findSmallDiagEntry(i,i+1);
- if (z==i-1)
+ Index j = findSmallDiagEntry(i,i+1);
+ if (j==i-1)
{
// block of (S T^{-1})
Matrix2s STi = m_T.template block<2,2>(i,i).template triangularView<Upper>().
@@ -352,7 +352,7 @@ namespace Eigen {
}
else
{
- pushDownZero(z,i,i+1);
+ pushDownZero(j,i,i+1);
}
}
diff --git a/Eigen/src/SVD/BDCSVD.h b/Eigen/src/SVD/BDCSVD.h
index 9b141c8df..55a83b574 100644
--- a/Eigen/src/SVD/BDCSVD.h
+++ b/Eigen/src/SVD/BDCSVD.h
@@ -825,7 +825,7 @@ void BDCSVD<MatrixType>::computeSingVals(const ArrayRef& col0, const ArrayRef& d
while (rightShifted - leftShifted > 2 * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(abs(leftShifted), abs(rightShifted)))
{
RealScalar midShifted = (leftShifted + rightShifted) / 2;
- RealScalar fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
+ fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
if (fLeft * fMid < 0)
{
rightShifted = midShifted;
diff --git a/Eigen/src/SVD/JacobiSVD.h b/Eigen/src/SVD/JacobiSVD.h
index a46a47104..e29d36cf2 100644
--- a/Eigen/src/SVD/JacobiSVD.h
+++ b/Eigen/src/SVD/JacobiSVD.h
@@ -387,7 +387,7 @@ struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
if(work_matrix.coeff(p,q) != Scalar(0))
{
- Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
+ z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
work_matrix.col(q) *= z;
if(svd.computeV()) svd.m_matrixV.col(q) *= z;
}
diff --git a/Eigen/src/SparseCore/SparseMatrix.h b/Eigen/src/SparseCore/SparseMatrix.h
index 5c5bf2d8c..d10880cb5 100644
--- a/Eigen/src/SparseCore/SparseMatrix.h
+++ b/Eigen/src/SparseCore/SparseMatrix.h
@@ -1119,9 +1119,9 @@ typename SparseMatrix<_Scalar,_Options,_Index>::Scalar& SparseMatrix<_Scalar,_Op
// so that the entire free-space is allocated to the current inner-vector.
eigen_internal_assert(data_end < m_data.allocatedSize());
StorageIndex new_end = convert_index(m_data.allocatedSize());
- for(Index j=outer+1; j<=m_outerSize; ++j)
- if(m_outerIndex[j]==data_end)
- m_outerIndex[j] = new_end;
+ for(Index k=outer+1; k<=m_outerSize; ++k)
+ if(m_outerIndex[k]==data_end)
+ m_outerIndex[k] = new_end;
}
return m_data.value(p);
}
@@ -1144,9 +1144,9 @@ typename SparseMatrix<_Scalar,_Options,_Index>::Scalar& SparseMatrix<_Scalar,_Op
// so that the entire free-space is allocated to the current inner-vector.
eigen_internal_assert(data_end < m_data.allocatedSize());
StorageIndex new_end = convert_index(m_data.allocatedSize());
- for(Index j=outer+1; j<=m_outerSize; ++j)
- if(m_outerIndex[j]==data_end)
- m_outerIndex[j] = new_end;
+ for(Index k=outer+1; k<=m_outerSize; ++k)
+ if(m_outerIndex[k]==data_end)
+ m_outerIndex[k] = new_end;
}
// and insert it at the right position (sorted insertion)
diff --git a/Eigen/src/SparseCore/SparseView.h b/Eigen/src/SparseCore/SparseView.h
index 1c69aa458..0a87f01d9 100644
--- a/Eigen/src/SparseCore/SparseView.h
+++ b/Eigen/src/SparseCore/SparseView.h
@@ -36,9 +36,9 @@ public:
EIGEN_SPARSE_PUBLIC_INTERFACE(SparseView)
typedef typename internal::remove_all<MatrixType>::type NestedExpression;
- explicit SparseView(const MatrixType& mat, const Scalar& m_reference = Scalar(0),
- RealScalar m_epsilon = NumTraits<Scalar>::dummy_precision()) :
- m_matrix(mat), m_reference(m_reference), m_epsilon(m_epsilon) {}
+ explicit SparseView(const MatrixType& mat, const Scalar& reference = Scalar(0),
+ RealScalar epsilon = NumTraits<Scalar>::dummy_precision())
+ : m_matrix(mat), m_reference(reference), m_epsilon(epsilon) {}
inline Index rows() const { return m_matrix.rows(); }
inline Index cols() const { return m_matrix.cols(); }
diff --git a/test/inverse.cpp b/test/inverse.cpp
index b09989aca..5c6777a18 100644
--- a/test/inverse.cpp
+++ b/test/inverse.cpp
@@ -71,11 +71,11 @@ template<typename MatrixType> void inverse(const MatrixType& m)
// check with submatrices
{
- Matrix<Scalar, MatrixType::RowsAtCompileTime+1, MatrixType::RowsAtCompileTime+1, MatrixType::Options> m3;
- m3.setRandom();
- m3.topLeftCorner(rows,rows) = m1;
- m2 = m3.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>().inverse();
- VERIFY_IS_APPROX( (m3.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>()), m2.inverse() );
+ Matrix<Scalar, MatrixType::RowsAtCompileTime+1, MatrixType::RowsAtCompileTime+1, MatrixType::Options> m5;
+ m5.setRandom();
+ m5.topLeftCorner(rows,rows) = m1;
+ m2 = m5.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>().inverse();
+ VERIFY_IS_APPROX( (m5.template topLeftCorner<MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime>()), m2.inverse() );
}
#endif
diff --git a/test/sparse_solver.h b/test/sparse_solver.h
index 418c18d6a..6ccc5487b 100644
--- a/test/sparse_solver.h
+++ b/test/sparse_solver.h
@@ -267,7 +267,6 @@ template<typename Solver> void check_sparse_spd_solving(Solver& solver)
{
if (it.sym() == SPD){
A = it.matrix();
- Mat halfA;
DenseVector b = it.rhs();
DenseVector refX = it.refX();
PermutationMatrix<Dynamic, Dynamic, StorageIndex> pnull;
diff --git a/test/sparse_vector.cpp b/test/sparse_vector.cpp
index 5dc421976..d3975b99b 100644
--- a/test/sparse_vector.cpp
+++ b/test/sparse_vector.cpp
@@ -55,16 +55,16 @@ template<typename Scalar,typename Index> void sparse_vector(int rows, int cols)
// test coeffRef with reallocation
{
- SparseVectorType v1(rows);
- DenseVector v2 = DenseVector::Zero(rows);
+ SparseVectorType v4(rows);
+ DenseVector v5 = DenseVector::Zero(rows);
for(int k=0; k<rows; ++k)
{
int i = internal::random<int>(0,rows-1);
Scalar v = internal::random<Scalar>();
- v1.coeffRef(i) += v;
- v2.coeffRef(i) += v;
+ v4.coeffRef(i) += v;
+ v5.coeffRef(i) += v;
}
- VERIFY_IS_APPROX(v1,v2);
+ VERIFY_IS_APPROX(v4,v5);
}
v1.coeffRef(nonzerocoords[0]) = Scalar(5);