Add cholesky's members to MatrixBase
Various documentation improvements including new snippets (AngleAxis and Cholesky)
diff --git a/Eigen/src/Cholesky/Cholesky.h b/Eigen/src/Cholesky/Cholesky.h
index 8d9b6e8..c1f05d7 100644
--- a/Eigen/src/Cholesky/Cholesky.h
+++ b/Eigen/src/Cholesky/Cholesky.h
@@ -66,7 +66,7 @@
     bool isPositiveDefinite(void) const { return m_isPositiveDefinite; }
 
     template<typename Derived>
-    typename Derived::Eval solve(MatrixBase<Derived> &b);
+    typename Derived::Eval solve(const MatrixBase<Derived> &b) const;
 
     void compute(const MatrixType& matrix);
 
@@ -110,10 +110,14 @@
 /** \returns the solution of A x = \a b using the current decomposition of A.
   * In other words, it returns \code A^-1 b \endcode computing
   * \code L^-*  L^1 b \endcode from right to left.
+  *
+  * Example: \include Cholesky_solve.cpp
+  * Output: \verbinclude Cholesky_solve.out
+  *
   */
 template<typename MatrixType>
 template<typename Derived>
-typename Derived::Eval Cholesky<MatrixType>::solve(MatrixBase<Derived> &b)
+typename Derived::Eval Cholesky<MatrixType>::solve(const MatrixBase<Derived> &b) const
 {
   const int size = m_matrix.rows();
   ei_assert(size==b.size());
@@ -121,5 +125,14 @@
   return m_matrix.adjoint().template extract<Upper>().inverseProduct(matrixL().inverseProduct(b));
 }
 
+/** \cholesky_module
+  * \returns the Cholesky decomposition of \c *this
+  */
+template<typename Derived>
+inline const Cholesky<typename ei_eval<Derived>::type>
+MatrixBase<Derived>::cholesky() const
+{
+  return Cholesky<typename ei_eval<Derived>::type>(derived());
+}
 
 #endif // EIGEN_CHOLESKY_H
diff --git a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h
index 8905385..2d85f78 100644
--- a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h
+++ b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h
@@ -77,7 +77,7 @@
     }
 
     template<typename Derived>
-    typename Derived::Eval solve(MatrixBase<Derived> &b);
+    typename Derived::Eval solve(const MatrixBase<Derived> &b) const;
 
     void compute(const MatrixType& matrix);
 
@@ -127,7 +127,7 @@
   */
 template<typename MatrixType>
 template<typename Derived>
-typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(MatrixBase<Derived> &vecB)
+typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(const MatrixBase<Derived> &vecB) const
 {
   const int size = m_matrix.rows();
   ei_assert(size==vecB.size());
@@ -140,5 +140,14 @@
       );
 }
 
+/** \cholesky_module
+  * \returns the Cholesky decomposition without square root of \c *this
+  */
+template<typename Derived>
+inline const CholeskyWithoutSquareRoot<typename ei_eval<Derived>::type>
+MatrixBase<Derived>::choleskyNoSqrt() const
+{
+  return derived();
+}
 
 #endif // EIGEN_CHOLESKY_WITHOUT_SQUARE_ROOT_H
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index c6ea5f1..bd4b64b 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -532,6 +532,10 @@
     void computeInverse(typename ei_eval<Derived>::type *result) const;
     Scalar determinant() const;
 
+/////////// Cholesky module ///////////
+
+    const Cholesky<typename ei_eval<Derived>::type> cholesky() const;
+    const CholeskyWithoutSquareRoot<typename ei_eval<Derived>::type> choleskyNoSqrt() const;
 
 /////////// QR module ///////////
 
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index 0246b43..1fba262 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -96,6 +96,8 @@
 
 template<typename ExpressionType, bool CheckExistence = true> class Inverse;
 template<typename MatrixType> class QR;
+template<typename MatrixType> class Cholesky;
+template<typename MatrixType> class CholeskyWithoutSquareRoot;
 
 // Geometry module:
 template<typename Lhs, typename Rhs> class Cross;
diff --git a/Eigen/src/Geometry/AngleAxis.h b/Eigen/src/Geometry/AngleAxis.h
index abb9b50..647e075 100644
--- a/Eigen/src/Geometry/AngleAxis.h
+++ b/Eigen/src/Geometry/AngleAxis.h
@@ -29,7 +29,7 @@
   *
   * \class AngleAxis
   *
-  * \brief Represents a 3D rotation as a rotation angle around an arbitray 3D axis
+  * \brief Represents a 3D rotation as a rotation angle around an arbitrary 3D axis
   *
   * \param _Scalar the scalar type, i.e., the type of the coefficients.
   *
@@ -37,7 +37,14 @@
   * \li \c AngleAxisf for \c float
   * \li \c AngleAxisd for \c double
   *
-  * \sa class Quaternion, class Transform
+  * \addexample AngleAxisForEuler \label How to define a rotation from Euler-angles
+  *
+  * Combined with MatrixBase::Unit{X,Y,Z}, AngleAxis can be used to easily
+  * mimic Euler-angles. Here is an example:
+  * \include AngleAxis_mimic_euler.cpp
+  * Output: \verbinclude AngleAxis_mimic_euler.out
+  *
+  * \sa class Quaternion, class Transform, MatrixBase::UnitX()
   */
 template<typename _Scalar>
 class AngleAxis
diff --git a/doc/Doxyfile.in b/doc/Doxyfile.in
index bade60a..4c5c71f 100644
--- a/doc/Doxyfile.in
+++ b/doc/Doxyfile.in
@@ -199,6 +199,7 @@
 ALIASES                = "only_for_vectors=This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column." \
                          "array_module=This is defined in the %Array module. \code #include <Eigen/Array> \endcode" \
                          "lu_module=This is defined in the %LU module. \code #include <Eigen/LU> \endcode" \
+                         "cholesky_module=This is defined in the %Cholesky module. \code #include <Eigen/Cholesky> \endcode" \
                          "qr_module=This is defined in the %QR module. \code #include <Eigen/QR> \endcode" \
                          "geometry_module=This is defined in the %Geometry module. \code #include <Eigen/Geometry> \endcode" \
                          "addexample=\anchor"  \
diff --git a/doc/buildexamplelist.sh b/doc/buildexamplelist.sh
index 0a9bb7e..7d47ae3 100755
--- a/doc/buildexamplelist.sh
+++ b/doc/buildexamplelist.sh
@@ -4,7 +4,7 @@
 echo "/** \page ExampleList"
 echo "<h1>Selected list of examples</h1>"
 
-grep \\addexample $1/Eigen/* -R | cut -d \\ -f 2- | \
+grep \\addexample $1/Eigen/src/*/*.h -R | cut -d \\ -f 2- | \
 while read example;
 do
 anchor=`echo "$example" | cut -d " " -f 2`
diff --git a/doc/snippets/AngleAxis_mimic_euler.cpp b/doc/snippets/AngleAxis_mimic_euler.cpp
new file mode 100644
index 0000000..be6b8ad
--- /dev/null
+++ b/doc/snippets/AngleAxis_mimic_euler.cpp
@@ -0,0 +1,4 @@
+Matrix3f m = AngleAxisf(0.25*M_PI, Vector3f::UnitX())
+           * AngleAxisf(0.5*M_PI,  Vector3f::UnitY())
+           * AngleAxisf(0.33*M_PI, Vector3f::UnitZ());
+cout << m << endl;
diff --git a/doc/snippets/CMakeLists.txt b/doc/snippets/CMakeLists.txt
index faf6440..72bd777 100644
--- a/doc/snippets/CMakeLists.txt
+++ b/doc/snippets/CMakeLists.txt
@@ -20,4 +20,6 @@
   ARGS >${CMAKE_CURRENT_BINARY_DIR}/${snippet}.out
 )
 ADD_DEPENDENCIES(all_snippets ${compile_snippet_target})
+set_source_files_properties(${CMAKE_CURRENT_BINARY_DIR}/${compile_snippet_src}
+                            PROPERTIES OBJECT_DEPENDS ${snippet_src})
 ENDFOREACH(snippet_src)
diff --git a/doc/snippets/Cholesky_solve.cpp b/doc/snippets/Cholesky_solve.cpp
new file mode 100644
index 0000000..4f9ac9c
--- /dev/null
+++ b/doc/snippets/Cholesky_solve.cpp
@@ -0,0 +1,6 @@
+typedef Matrix<float,Dynamic,2> DataMatrix;
+// let's generate some samples on the 3D plane of equation z = 2x+3y (with some noise)
+DataMatrix samples = DataMatrix::random(12,2);
+VectorXf elevations = 2*samples.col(0) + 3*samples.col(1) + VectorXf::random(12)*0.1;
+// and let's solve samples * x = elevations in least square sense:
+cout << (samples.adjoint() * samples).cholesky().solve((samples.adjoint()*elevations).eval()) << endl;
diff --git a/doc/snippets/compile_snippet.cpp.in b/doc/snippets/compile_snippet.cpp.in
index 5876aab..950f066 100644
--- a/doc/snippets/compile_snippet.cpp.in
+++ b/doc/snippets/compile_snippet.cpp.in
@@ -1,8 +1,11 @@
 #include <Eigen/Core>
 #include <Eigen/Array>
 #include <Eigen/LU>
+#include <Eigen/Cholesky>
+#include <Eigen/Geometry>
 
 USING_PART_OF_NAMESPACE_EIGEN
+using namespace Eigen;
 using namespace std;
 
 int main(int, char**)