aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/CMakeLists.txt2
-rw-r--r--Eigen/Cholesky2
-rw-r--r--Eigen/Regression5
-rw-r--r--Eigen/src/Cholesky/Cholesky.h165
-rw-r--r--Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h184
-rw-r--r--Eigen/src/Core/DiagonalMatrix.h7
-rw-r--r--Eigen/src/Core/Dot.h16
-rw-r--r--Eigen/src/Core/MatrixBase.h13
-rw-r--r--Eigen/src/Core/Minor.h9
-rw-r--r--Eigen/src/Core/Part.h9
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h3
-rw-r--r--Eigen/src/Core/util/Memory.h6
-rw-r--r--Eigen/src/Geometry/Transform.h10
-rw-r--r--disabled/Householder.h5
-rw-r--r--doc/Doxyfile.in2
-rw-r--r--doc/Experimental.dox47
16 files changed, 81 insertions, 404 deletions
diff --git a/Eigen/CMakeLists.txt b/Eigen/CMakeLists.txt
index 7eeb595dd..76aca5e4d 100644
--- a/Eigen/CMakeLists.txt
+++ b/Eigen/CMakeLists.txt
@@ -1,4 +1,4 @@
-set(Eigen_HEADERS Core LU Cholesky QR Geometry Sparse Array SVD Regression LeastSquares StdVector)
+set(Eigen_HEADERS Core LU Cholesky QR Geometry Sparse Array SVD LeastSquares StdVector)
if(EIGEN_BUILD_LIB)
set(Eigen_SRCS
diff --git a/Eigen/Cholesky b/Eigen/Cholesky
index 3c3ffb936..b2179fff9 100644
--- a/Eigen/Cholesky
+++ b/Eigen/Cholesky
@@ -31,8 +31,6 @@ namespace Eigen {
#include "src/Array/Functors.h"
#include "src/Cholesky/LLT.h"
#include "src/Cholesky/LDLT.h"
-#include "src/Cholesky/Cholesky.h"
-#include "src/Cholesky/CholeskyWithoutSquareRoot.h"
} // namespace Eigen
diff --git a/Eigen/Regression b/Eigen/Regression
deleted file mode 100644
index cc6bae89b..000000000
--- a/Eigen/Regression
+++ /dev/null
@@ -1,5 +0,0 @@
-#ifdef __GNUC__
-#warning "The Eigen/Regression header file has been renamed to Eigen/LeastSquares. The old name is deprecated, please update your code."
-#endif
-
-#include "LeastSquares" \ No newline at end of file
diff --git a/Eigen/src/Cholesky/Cholesky.h b/Eigen/src/Cholesky/Cholesky.h
deleted file mode 100644
index 31451a6b0..000000000
--- a/Eigen/src/Cholesky/Cholesky.h
+++ /dev/null
@@ -1,165 +0,0 @@
-// This file is part of Eigen, a lightweight C++ template library
-// for linear algebra. Eigen itself is part of the KDE project.
-//
-// Copyright (C) 2008 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_CHOLESKY_H
-#define EIGEN_CHOLESKY_H
-
-/** \ingroup Cholesky_Module
- *
- * \class Cholesky
- *
- * \deprecated this class has been renamed LLT
- */
-template<typename MatrixType> class Cholesky
-{
- private:
- typedef typename MatrixType::Scalar Scalar;
- typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
-
- enum {
- PacketSize = ei_packet_traits<Scalar>::size,
- AlignmentMask = int(PacketSize)-1
- };
-
- public:
-
- Cholesky(const MatrixType& matrix)
- : m_matrix(matrix.rows(), matrix.cols())
- {
- compute(matrix);
- }
-
- /** \deprecated */
- inline Part<MatrixType, LowerTriangular> matrixL(void) const { return m_matrix; }
-
- /** \deprecated */
- inline bool isPositiveDefinite(void) const { return m_isPositiveDefinite; }
-
- template<typename Derived>
- EIGEN_DEPRECATED typename MatrixBase<Derived>::PlainMatrixType_ColMajor solve(const MatrixBase<Derived> &b) const;
-
- template<typename RhsDerived, typename ResDerived>
- bool solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const;
-
- template<typename Derived>
- bool solveInPlace(MatrixBase<Derived> &bAndX) const;
-
- void compute(const MatrixType& matrix);
-
- protected:
- /** \internal
- * Used to compute and store L
- * The strict upper part is not used and even not initialized.
- */
- MatrixType m_matrix;
- bool m_isPositiveDefinite;
-};
-
-/** \deprecated */
-template<typename MatrixType>
-void Cholesky<MatrixType>::compute(const MatrixType& a)
-{
- assert(a.rows()==a.cols());
- const int size = a.rows();
- m_matrix.resize(size, size);
- const RealScalar eps = ei_sqrt(precision<Scalar>());
-
- RealScalar x;
- x = ei_real(a.coeff(0,0));
- m_isPositiveDefinite = x > eps && ei_isMuchSmallerThan(ei_imag(a.coeff(0,0)), RealScalar(1));
- m_matrix.coeffRef(0,0) = ei_sqrt(x);
- m_matrix.col(0).end(size-1) = a.row(0).end(size-1).adjoint() / ei_real(m_matrix.coeff(0,0));
- for (int j = 1; j < size; ++j)
- {
- Scalar tmp = ei_real(a.coeff(j,j)) - m_matrix.row(j).start(j).squaredNorm();
- x = ei_real(tmp);
- if (x < eps || (!ei_isMuchSmallerThan(ei_imag(tmp), RealScalar(1))))
- {
- m_isPositiveDefinite = false;
- return;
- }
- m_matrix.coeffRef(j,j) = x = ei_sqrt(x);
-
- int endSize = size-j-1;
- if (endSize>0) {
- // Note that when all matrix columns have good alignment, then the following
- // product is guaranteed to be optimal with respect to alignment.
- m_matrix.col(j).end(endSize) =
- (m_matrix.block(j+1, 0, endSize, j) * m_matrix.row(j).start(j).adjoint()).lazy();
-
- // FIXME could use a.col instead of a.row
- m_matrix.col(j).end(endSize) = (a.row(j).end(endSize).adjoint()
- - m_matrix.col(j).end(endSize) ) / x;
- }
- }
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename Derived>
-typename MatrixBase<Derived>::PlainMatrixType_ColMajor Cholesky<MatrixType>::solve(const MatrixBase<Derived> &b) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==b.rows());
- typename MatrixBase<Derived>::PlainMatrixType_ColMajor x(b);
- solveInPlace(x);
- return x;
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename RhsDerived, typename ResDerived>
-bool Cholesky<MatrixType>::solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==b.rows() && "Cholesky::solve(): invalid number of rows of the right hand side matrix b");
- return solveInPlace((*result) = b);
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename Derived>
-bool Cholesky<MatrixType>::solveInPlace(MatrixBase<Derived> &bAndX) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==bAndX.rows());
- if (!m_isPositiveDefinite)
- return false;
- matrixL().solveTriangularInPlace(bAndX);
- m_matrix.adjoint().template part<UpperTriangular>().solveTriangularInPlace(bAndX);
- return true;
-}
-
-/** \cholesky_module
- * \deprecated has been renamed llt()
- */
-template<typename Derived>
-inline const Cholesky<typename MatrixBase<Derived>::PlainMatrixType>
-MatrixBase<Derived>::cholesky() const
-{
- return Cholesky<PlainMatrixType>(derived());
-}
-
-#endif // EIGEN_CHOLESKY_H
diff --git a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h
deleted file mode 100644
index 641f9be01..000000000
--- a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h
+++ /dev/null
@@ -1,184 +0,0 @@
-// This file is part of Eigen, a lightweight C++ template library
-// for linear algebra. Eigen itself is part of the KDE project.
-//
-// Copyright (C) 2008 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_CHOLESKY_WITHOUT_SQUARE_ROOT_H
-#define EIGEN_CHOLESKY_WITHOUT_SQUARE_ROOT_H
-
-/** \deprecated \ingroup Cholesky_Module
- *
- * \class CholeskyWithoutSquareRoot
- *
- * \deprecated this class has been renamed LDLT
- */
-template<typename MatrixType> class CholeskyWithoutSquareRoot
-{
- public:
-
- typedef typename MatrixType::Scalar Scalar;
- typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
-
- CholeskyWithoutSquareRoot(const MatrixType& matrix)
- : m_matrix(matrix.rows(), matrix.cols())
- {
- compute(matrix);
- }
-
- /** \returns the lower triangular matrix L */
- inline Part<MatrixType, UnitLowerTriangular> matrixL(void) const { return m_matrix; }
-
- /** \returns the coefficients of the diagonal matrix D */
- inline DiagonalCoeffs<MatrixType> vectorD(void) const { return m_matrix.diagonal(); }
-
- /** \returns true if the matrix is positive definite */
- inline bool isPositiveDefinite(void) const { return m_isPositiveDefinite; }
-
- template<typename Derived>
- EIGEN_DEPRECATED typename Derived::Eval solve(const MatrixBase<Derived> &b) const;
-
- template<typename RhsDerived, typename ResDerived>
- bool solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const;
-
- template<typename Derived>
- bool solveInPlace(MatrixBase<Derived> &bAndX) const;
-
- void compute(const MatrixType& matrix);
-
- protected:
- /** \internal
- * Used to compute and store the cholesky decomposition A = L D L^* = U^* D U.
- * The strict upper part is used during the decomposition, the strict lower
- * part correspond to the coefficients of L (its diagonal is equal to 1 and
- * is not stored), and the diagonal entries correspond to D.
- */
- MatrixType m_matrix;
-
- bool m_isPositiveDefinite;
-};
-
-/** \deprecated */
-template<typename MatrixType>
-void CholeskyWithoutSquareRoot<MatrixType>::compute(const MatrixType& a)
-{
- assert(a.rows()==a.cols());
- const int size = a.rows();
- m_matrix.resize(size, size);
- m_isPositiveDefinite = true;
- const RealScalar eps = ei_sqrt(precision<Scalar>());
-
- if (size<=1)
- {
- m_matrix = a;
- return;
- }
-
- // Let's preallocate a temporay vector to evaluate the matrix-vector product into it.
- // Unlike the standard Cholesky decomposition, here we cannot evaluate it to the destination
- // matrix because it a sub-row which is not compatible suitable for efficient packet evaluation.
- // (at least if we assume the matrix is col-major)
- Matrix<Scalar,MatrixType::RowsAtCompileTime,1> _temporary(size);
-
- // Note that, in this algorithm the rows of the strict upper part of m_matrix is used to store
- // column vector, thus the strange .conjugate() and .transpose()...
-
- m_matrix.row(0) = a.row(0).conjugate();
- m_matrix.col(0).end(size-1) = m_matrix.row(0).end(size-1) / m_matrix.coeff(0,0);
- for (int j = 1; j < size; ++j)
- {
- RealScalar tmp = ei_real(a.coeff(j,j) - (m_matrix.row(j).start(j) * m_matrix.col(j).start(j).conjugate()).coeff(0,0));
- m_matrix.coeffRef(j,j) = tmp;
-
- if (tmp < eps)
- {
- m_isPositiveDefinite = false;
- return;
- }
-
- int endSize = size-j-1;
- if (endSize>0)
- {
- _temporary.end(endSize) = ( m_matrix.block(j+1,0, endSize, j)
- * m_matrix.col(j).start(j).conjugate() ).lazy();
-
- m_matrix.row(j).end(endSize) = a.row(j).end(endSize).conjugate()
- - _temporary.end(endSize).transpose();
-
- m_matrix.col(j).end(endSize) = m_matrix.row(j).end(endSize) / tmp;
- }
- }
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename Derived>
-typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(const MatrixBase<Derived> &b) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==b.rows());
-
- return m_matrix.adjoint().template part<UnitUpperTriangular>()
- .solveTriangular(
- ( m_matrix.cwise().inverse().template part<Diagonal>()
- * matrixL().solveTriangular(b))
- );
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename RhsDerived, typename ResDerived>
-bool CholeskyWithoutSquareRoot<MatrixType>
-::solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==b.rows() && "Cholesky::solve(): invalid number of rows of the right hand side matrix b");
- *result = b;
- return solveInPlace(*result);
-}
-
-/** \deprecated */
-template<typename MatrixType>
-template<typename Derived>
-bool CholeskyWithoutSquareRoot<MatrixType>::solveInPlace(MatrixBase<Derived> &bAndX) const
-{
- const int size = m_matrix.rows();
- ei_assert(size==bAndX.rows());
- if (!m_isPositiveDefinite)
- return false;
- matrixL().solveTriangularInPlace(bAndX);
- bAndX = (m_matrix.cwise().inverse().template part<Diagonal>() * bAndX).lazy();
- m_matrix.adjoint().template part<UnitUpperTriangular>().solveTriangularInPlace(bAndX);
- return true;
-}
-
-/** \cholesky_module
- * \deprecated has been renamed ldlt()
- */
-template<typename Derived>
-inline const CholeskyWithoutSquareRoot<typename MatrixBase<Derived>::PlainMatrixType>
-MatrixBase<Derived>::choleskyNoSqrt() const
-{
- return derived();
-}
-
-#endif // EIGEN_CHOLESKY_WITHOUT_SQUARE_ROOT_H
diff --git a/Eigen/src/Core/DiagonalMatrix.h b/Eigen/src/Core/DiagonalMatrix.h
index 25fa26953..07eaf0747 100644
--- a/Eigen/src/Core/DiagonalMatrix.h
+++ b/Eigen/src/Core/DiagonalMatrix.h
@@ -26,6 +26,7 @@
#define EIGEN_DIAGONALMATRIX_H
/** \class DiagonalMatrix
+ * \nonstableyet
*
* \brief Expression of a diagonal matrix
*
@@ -91,7 +92,8 @@ class DiagonalMatrix : ei_no_assignment_operator,
const typename CoeffsVectorType::Nested m_coeffs;
};
-/** \returns an expression of a diagonal matrix with *this as vector of diagonal coefficients
+/** \nonstableyet
+ * \returns an expression of a diagonal matrix with *this as vector of diagonal coefficients
*
* \only_for_vectors
*
@@ -109,7 +111,8 @@ MatrixBase<Derived>::asDiagonal() const
return derived();
}
-/** \returns true if *this is approximately equal to a diagonal matrix,
+/** \nonstableyet
+ * \returns true if *this is approximately equal to a diagonal matrix,
* within the precision given by \a prec.
*
* Example: \include MatrixBase_isDiagonal.cpp
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index b66fcbaae..4263a27c5 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -271,22 +271,6 @@ MatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
/** \returns the squared norm of *this, i.e. the dot product of *this with itself.
*
- * \note This is \em not the \em l2 norm, but its square.
- *
- * \deprecated Use squaredNorm() instead. This norm2() function is kept only for compatibility and will be removed in Eigen 2.0.
- *
- * \only_for_vectors
- *
- * \sa dot(), norm()
- */
-template<typename Derived>
-EIGEN_DEPRECATED inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real MatrixBase<Derived>::norm2() const
-{
- return ei_real((*this).cwise().abs2().sum());
-}
-
-/** \returns the squared norm of *this, i.e. the dot product of *this with itself.
- *
* \only_for_vectors
*
* \sa dot(), norm()
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index d342a8936..8c755c592 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -229,10 +229,6 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
Derived& operator=(const MatrixBase<OtherDerived>& other);
- /** Copies \a other into *this without evaluating other. \returns a reference to *this. */
- template<typename OtherDerived>
- Derived& lazyAssign(const MatrixBase<OtherDerived>& other);
-
/** Special case of the template operator=, in order to prevent the compiler
* from generating a default operator= (issue hit with g++ 4.1)
*/
@@ -241,6 +237,11 @@ template<typename Derived> class MatrixBase
return this->operator=<Derived>(other);
}
+#ifndef EIGEN_PARSED_BY_DOXYGEN
+ /** Copies \a other into *this without evaluating other. \returns a reference to *this. */
+ template<typename OtherDerived>
+ Derived& lazyAssign(const MatrixBase<OtherDerived>& other);
+
/** Overloaded for cache friendly product evaluation */
template<typename Lhs, typename Rhs>
Derived& lazyAssign(const Product<Lhs,Rhs,CacheFriendlyProduct>& product);
@@ -249,6 +250,7 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
Derived& lazyAssign(const Flagged<OtherDerived, 0, EvalBeforeNestingBit | EvalBeforeAssigningBit>& other)
{ return lazyAssign(other._expression()); }
+#endif // not EIGEN_PARSED_BY_DOXYGEN
CommaInitializer<Derived> operator<< (const Scalar& s);
@@ -589,9 +591,6 @@ template<typename Derived> class MatrixBase
const LLT<PlainMatrixType> llt() const;
const LDLT<PlainMatrixType> ldlt() const;
- // deprecated:
- const Cholesky<PlainMatrixType> cholesky() const;
- const CholeskyWithoutSquareRoot<PlainMatrixType> choleskyNoSqrt() const;
/////////// QR module ///////////
diff --git a/Eigen/src/Core/Minor.h b/Eigen/src/Core/Minor.h
index 829bfefb2..e2d47da79 100644
--- a/Eigen/src/Core/Minor.h
+++ b/Eigen/src/Core/Minor.h
@@ -25,7 +25,8 @@
#ifndef EIGEN_MINOR_H
#define EIGEN_MINOR_H
-/** \class Minor
+/** \nonstableyet
+ * \class Minor
*
* \brief Expression of a minor
*
@@ -92,7 +93,8 @@ template<typename MatrixType> class Minor
const int m_row, m_col;
};
-/** \return an expression of the (\a row, \a col)-minor of *this,
+/** \nonstableyet
+ * \return an expression of the (\a row, \a col)-minor of *this,
* i.e. an expression constructed from *this by removing the specified
* row and column.
*
@@ -108,7 +110,8 @@ MatrixBase<Derived>::minor(int row, int col)
return Minor<Derived>(derived(), row, col);
}
-/** This is the const version of minor(). */
+/** \nonstableyet
+ * This is the const version of minor(). */
template<typename Derived>
inline const Minor<Derived>
MatrixBase<Derived>::minor(int row, int col) const
diff --git a/Eigen/src/Core/Part.h b/Eigen/src/Core/Part.h
index b0ef1ff3d..d79c2f2ad 100644
--- a/Eigen/src/Core/Part.h
+++ b/Eigen/src/Core/Part.h
@@ -26,7 +26,8 @@
#ifndef EIGEN_PART_H
#define EIGEN_PART_H
-/** \class Part
+/** \nonstableyet
+ * \class Part
*
* \brief Expression of a triangular matrix extracted from a given matrix
*
@@ -127,7 +128,8 @@ template<typename MatrixType, unsigned int Mode> class Part
const typename MatrixType::Nested m_matrix;
};
-/** \returns an expression of a triangular matrix extracted from the current matrix
+/** \nonstableyet
+ * \returns an expression of a triangular matrix extracted from the current matrix
*
* The parameter \a Mode can have the following values: \c UpperTriangular, \c StrictlyUpperTriangular, \c UnitUpperTriangular,
* \c LowerTriangular, \c StrictlyLowerTriangular, \c UnitLowerTriangular.
@@ -278,7 +280,8 @@ void Part<MatrixType, Mode>::lazyAssign(const Other& other)
>::run(m_matrix.const_cast_derived(), other.derived());
}
-/** \returns a lvalue pseudo-expression allowing to perform special operations on \c *this.
+/** \nonstableyet
+ * \returns a lvalue pseudo-expression allowing to perform special operations on \c *this.
*
* The \a Mode parameter can have the following values: \c UpperTriangular, \c StrictlyUpperTriangular, \c LowerTriangular,
* \c StrictlyLowerTriangular, \c SelfAdjoint.
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index a45210e0c..a72a40b1b 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -106,9 +106,6 @@ template<typename MatrixType> class QR;
template<typename MatrixType> class SVD;
template<typename MatrixType> class LLT;
template<typename MatrixType> class LDLT;
-// deprecated:
-template<typename MatrixType> class Cholesky;
-template<typename MatrixType> class CholeskyWithoutSquareRoot;
// Geometry module:
template<typename Derived, int _Dim> class RotationBase;
diff --git a/Eigen/src/Core/util/Memory.h b/Eigen/src/Core/util/Memory.h
index ae529797c..3c9fdf5cc 100644
--- a/Eigen/src/Core/util/Memory.h
+++ b/Eigen/src/Core/util/Memory.h
@@ -307,12 +307,6 @@ inline static int ei_alignmentOffset(const Scalar* ptr, int maxOffset)
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(Scalar,Size) \
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(((Size)!=Eigen::Dynamic) && ((sizeof(Scalar)*(Size))%16==0))
-/** Deprecated, use the EIGEN_MAKE_ALIGNED_OPERATOR_NEW macro instead in your own class */
-struct WithAlignedOperatorNew
-{
- EIGEN_MAKE_ALIGNED_OPERATOR_NEW
-};
-
/** \class aligned_allocator
*
* \brief stl compatible allocator to use with with 16 byte aligned types
diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h
index 83bad3afd..a932b3ff8 100644
--- a/Eigen/src/Geometry/Transform.h
+++ b/Eigen/src/Geometry/Transform.h
@@ -247,7 +247,6 @@ public:
template<typename Derived>
inline Transform operator*(const RotationBase<Derived,Dim>& r) const;
- EIGEN_DEPRECATED LinearMatrixType extractRotation(TransformTraits traits = Affine) const { return rotation(traits); }
LinearMatrixType rotation(TransformTraits traits = Affine) const;
template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
@@ -595,6 +594,7 @@ inline Transform<Scalar,Dim> Transform<Scalar,Dim>::operator*(const RotationBase
***************************/
/** \returns the rotation part of the transformation
+ * \nonstableyet
*
* \param traits allows to optimize the extraction process when the transformion
* is known to be not a general aafine transformation. The possible values are:
@@ -625,10 +625,10 @@ Transform<Scalar,Dim>::rotation(TransformTraits traits) const
return matQ;
}
else if (traits == Isometry) // though that's stupid let's handle it !
- return linear();
+ return linear(); // FIXME needs to divide by determinant
else
{
- ei_assert("invalid traits value in Transform::extractRotation()");
+ ei_assert("invalid traits value in Transform::rotation()");
return LinearMatrixType();
}
}
@@ -650,7 +650,9 @@ Transform<Scalar,Dim>::fromPositionOrientationScale(const MatrixBase<PositionDer
return *this;
}
-/** \returns the inverse transformation matrix according to some given knowledge
+/** \nonstableyet
+ *
+ * \returns the inverse transformation matrix according to some given knowledge
* on \c *this.
*
* \param traits allows to optimize the inversion process when the transformion
diff --git a/disabled/Householder.h b/disabled/Householder.h
index 43945b01a..874b812db 100644
--- a/disabled/Householder.h
+++ b/disabled/Householder.h
@@ -8,10 +8,11 @@ template<typename InputVector, typename OutputVector>
void ei_compute_householder(const InputVector& x, OutputVector *v, typename OutputVector::RealScalar *beta)
{
EIGEN_STATIC_ASSERT(ei_is_same_type<typename InputVector::Scalar, typename OutputVector::Scalar>::ret,
- YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
+ YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
EIGEN_STATIC_ASSERT((InputVector::SizeAtCompileTime == OutputVector::SizeAtCompileTime+1)
|| InputVector::SizeAtCompileTime == Dynamic
- || OutputVector::SizeAtCompileTime == Dynamic)
+ || OutputVector::SizeAtCompileTime == Dynamic,
+ YOU_MIXED_VECTORS_OF_DIFFERENT_SIZES)
typedef typename OutputVector::RealScalar RealScalar;
ei_assert(x.size() == v->size()+1);
int n = x.size();
diff --git a/doc/Doxyfile.in b/doc/Doxyfile.in
index 46187f4a2..b30ffe0d8 100644
--- a/doc/Doxyfile.in
+++ b/doc/Doxyfile.in
@@ -214,7 +214,7 @@ ALIASES = "only_for_vectors=This is only for vectors (either row-
"addexample=\anchor" \
"label=\bug" \
"redstar=<a href='#warningarraymodule' style='color:red;text-decoration: none;'>*</a>" \
- "nonstableyet=\warning This class/function is not considered to be part of the stable public API yet. Some (minor) changes might happen in future releases."
+ "nonstableyet=\warning This class/function is not considered to be part of the stable public API yet. Changes may happen in future releases. See \ref Experimental \"Experimental parts of Eigen\""
# Set the OPTIMIZE_OUTPUT_FOR_C tag to YES if your project consists of C
# sources only. Doxygen will then generate output that is more tailored for C.
diff --git a/doc/Experimental.dox b/doc/Experimental.dox
new file mode 100644
index 000000000..f55a47f50
--- /dev/null
+++ b/doc/Experimental.dox
@@ -0,0 +1,47 @@
+namespace Eigen {
+
+/** \page Experimental Experimental parts of Eigen
+
+\b Table \b of \b contents
+ - \ref summary
+ - \ref modules
+ - \ref core
+
+\section summary Summary
+
+Experimental features may at any time:
+\li be removed;
+\li be subject to an API incompatible change;
+\li introduce API or ABI incompatible changes in your own application if you let them affect your API or ABI.
+
+\section modules Experimental modules
+
+The following modules are considered entirely experimental:
+\li SVD
+\li QR
+\li Sparse
+
+\section core Experimental parts of the Core module
+
+In the Core module, the only classes subject to ABI stability guarantee (meaning that you can use it for data members in your public ABI) is:
+\li Matrix
+\li Map
+
+All other classes offer no ABI guarantee, e.g. the layout of their data can be changed.
+
+The only classes subject to (even partial) API stability guarantee (meaning that you can safely construct and use objects) are:
+\li MatrixBase : partial API stability (see below)
+\li Matrix : full API stability (except for experimental stuff inherited from MatrixBase)
+\li Map : full API stability (except for experimental stuff inherited from MatrixBase)
+
+All other classes offer no direct API guarantee, e.g. their methods can be changed; however notice that most classes inherit MatrixBase and that this is where most of their API comes from -- so in practice most of the API is stable.
+
+Here are the MatrixBase methods that are considered experimental, hence not part of any API stability guarantee:
+\li all methods documented as internal
+\li all methods hidden in the Doxygen documentation
+\li all methods marked as experimental
+\li all methods defined in experimental modules
+
+*/
+
+}