aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/QR
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-04-22 14:11:18 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-04-22 14:11:18 -0400
commit9962c59b56960569c8df332144190e62c1eb3b01 (patch)
treea3efa574460c6a08f4ed17a3896b497d5bfc374f /Eigen/src/QR
parent28dde19e40a3d758faa94f0fc228857f7b3192ea (diff)
* implement the corner() API change: new methods topLeftCorner() etc
* get rid of BlockReturnType: it was not needed, and code was not always using it consistently anyway * add topRows(), leftCols(), bottomRows(), rightCols() * add corners unit-test covering all of that * adapt docs, expand "porting from eigen 2 to 3" * adapt Eigen2Support
Diffstat (limited to 'Eigen/src/QR')
-rw-r--r--Eigen/src/QR/ColPivHouseholderQR.h14
-rw-r--r--Eigen/src/QR/FullPivHouseholderQR.h14
-rw-r--r--Eigen/src/QR/HouseholderQR.h12
3 files changed, 20 insertions, 20 deletions
diff --git a/Eigen/src/QR/ColPivHouseholderQR.h b/Eigen/src/QR/ColPivHouseholderQR.h
index 5893125cd..8c9c8840b 100644
--- a/Eigen/src/QR/ColPivHouseholderQR.h
+++ b/Eigen/src/QR/ColPivHouseholderQR.h
@@ -411,7 +411,7 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
{
m_nonzero_pivots = k;
m_hCoeffs.tail(size-k).setZero();
- m_qr.corner(BottomRight,rows-k,cols-k)
+ m_qr.bottomRightCorner(rows-k,cols-k)
.template triangularView<StrictlyLower>()
.setZero();
break;
@@ -436,7 +436,7 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
if(ei_abs(beta) > m_maxpivot) m_maxpivot = ei_abs(beta);
// apply the householder transformation
- m_qr.corner(BottomRight, rows-k, cols-k-1)
+ m_qr.bottomRightCorner(rows-k, cols-k-1)
.applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
// update our table of squared norms of the columns
@@ -483,17 +483,17 @@ struct ei_solve_retval<ColPivHouseholderQR<_MatrixType>, Rhs>
));
dec().matrixQR()
- .corner(TopLeft, nonzero_pivots, nonzero_pivots)
+ .topLeftCorner(nonzero_pivots, nonzero_pivots)
.template triangularView<Upper>()
- .solveInPlace(c.corner(TopLeft, nonzero_pivots, c.cols()));
+ .solveInPlace(c.topRows(nonzero_pivots));
typename Rhs::PlainObject d(c);
- d.corner(TopLeft, nonzero_pivots, c.cols())
+ d.topRows(nonzero_pivots)
= dec().matrixQR()
- .corner(TopLeft, nonzero_pivots, nonzero_pivots)
+ .topLeftCorner(nonzero_pivots, nonzero_pivots)
.template triangularView<Upper>()
- * c.corner(TopLeft, nonzero_pivots, c.cols());
+ * c.topRows(nonzero_pivots);
for(int i = 0; i < nonzero_pivots; ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
for(int i = nonzero_pivots; i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
diff --git a/Eigen/src/QR/FullPivHouseholderQR.h b/Eigen/src/QR/FullPivHouseholderQR.h
index 8df2ed1e0..0195e1330 100644
--- a/Eigen/src/QR/FullPivHouseholderQR.h
+++ b/Eigen/src/QR/FullPivHouseholderQR.h
@@ -314,7 +314,7 @@ FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(cons
int row_of_biggest_in_corner, col_of_biggest_in_corner;
RealScalar biggest_in_corner;
- biggest_in_corner = m_qr.corner(Eigen::BottomRight, rows-k, cols-k)
+ biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
.cwiseAbs()
.maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
row_of_biggest_in_corner += k;
@@ -349,7 +349,7 @@ FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(cons
m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
m_qr.coeffRef(k,k) = beta;
- m_qr.corner(BottomRight, rows-k, cols-k-1)
+ m_qr.bottomRightCorner(rows-k, cols-k-1)
.applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
}
@@ -389,7 +389,7 @@ struct ei_solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
{
int remainingSize = rows-k;
c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k)));
- c.corner(BottomRight, remainingSize, rhs().cols())
+ c.bottomRightCorner(remainingSize, rhs().cols())
.applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1),
dec().hCoeffs().coeff(k), &temp.coeffRef(0));
}
@@ -397,17 +397,17 @@ struct ei_solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
if(!dec().isSurjective())
{
// is c is in the image of R ?
- RealScalar biggest_in_upper_part_of_c = c.corner(TopLeft, dec().rank(), c.cols()).cwiseAbs().maxCoeff();
- RealScalar biggest_in_lower_part_of_c = c.corner(BottomLeft, rows-dec().rank(), c.cols()).cwiseAbs().maxCoeff();
+ RealScalar biggest_in_upper_part_of_c = c.topRows( dec().rank() ).cwiseAbs().maxCoeff();
+ RealScalar biggest_in_lower_part_of_c = c.bottomRows(rows-dec().rank()).cwiseAbs().maxCoeff();
// FIXME brain dead
const RealScalar m_precision = NumTraits<Scalar>::epsilon() * std::min(rows,cols);
if(!ei_isMuchSmallerThan(biggest_in_lower_part_of_c, biggest_in_upper_part_of_c, m_precision))
return;
}
dec().matrixQR()
- .corner(TopLeft, dec().rank(), dec().rank())
+ .topLeftCorner(dec().rank(), dec().rank())
.template triangularView<Upper>()
- .solveInPlace(c.corner(TopLeft, dec().rank(), c.cols()));
+ .solveInPlace(c.topRows(dec().rank()));
for(int i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
for(int i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
diff --git a/Eigen/src/QR/HouseholderQR.h b/Eigen/src/QR/HouseholderQR.h
index 10b9fb93f..6a2883939 100644
--- a/Eigen/src/QR/HouseholderQR.h
+++ b/Eigen/src/QR/HouseholderQR.h
@@ -216,7 +216,7 @@ HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const MatrixType&
m_qr.coeffRef(k,k) = beta;
// apply H to remaining part of m_qr from the left
- m_qr.corner(BottomRight, remainingRows, remainingCols)
+ m_qr.bottomRightCorner(remainingRows, remainingCols)
.applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingRows-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
}
m_isInitialized = true;
@@ -239,17 +239,17 @@ struct ei_solve_retval<HouseholderQR<_MatrixType>, Rhs>
// Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
c.applyOnTheLeft(householderSequence(
- dec().matrixQR().corner(TopLeft,rows,rank),
+ dec().matrixQR().leftCols(rank),
dec().hCoeffs().head(rank)).transpose()
);
dec().matrixQR()
- .corner(TopLeft, rank, rank)
+ .topLeftCorner(rank, rank)
.template triangularView<Upper>()
- .solveInPlace(c.corner(TopLeft, rank, c.cols()));
+ .solveInPlace(c.topRows(rank));
- dst.corner(TopLeft, rank, c.cols()) = c.corner(TopLeft, rank, c.cols());
- dst.corner(BottomLeft, cols-rank, c.cols()).setZero();
+ dst.topRows(rank) = c.topRows(rank);
+ dst.bottomRows(cols-rank).setZero();
}
};