aboutsummaryrefslogtreecommitdiffhomepage
path: root/doc/snippets
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2008-07-19 22:59:05 +0000
committerGravatar Gael Guennebaud <g.gael@free.fr>2008-07-19 22:59:05 +0000
commit269f683902231020a910fd4e2c1b74554183e2c8 (patch)
treef95a531cd1792f1df2122a3879217adc8eb9ec18 /doc/snippets
parent6e2c53e0562785f2f8c72a67dd8c40451cf2a319 (diff)
Add cholesky's members to MatrixBase
Various documentation improvements including new snippets (AngleAxis and Cholesky)
Diffstat (limited to 'doc/snippets')
-rw-r--r--doc/snippets/AngleAxis_mimic_euler.cpp4
-rw-r--r--doc/snippets/CMakeLists.txt2
-rw-r--r--doc/snippets/Cholesky_solve.cpp6
-rw-r--r--doc/snippets/compile_snippet.cpp.in3
4 files changed, 15 insertions, 0 deletions
diff --git a/doc/snippets/AngleAxis_mimic_euler.cpp b/doc/snippets/AngleAxis_mimic_euler.cpp
new file mode 100644
index 000000000..be6b8adbe
--- /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 faf6440e8..72bd7770a 100644
--- a/doc/snippets/CMakeLists.txt
+++ b/doc/snippets/CMakeLists.txt
@@ -20,4 +20,6 @@ ADD_CUSTOM_COMMAND(
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 000000000..4f9ac9c1e
--- /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 5876aab9c..950f06654 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**)