aboutsummaryrefslogtreecommitdiffhomepage
path: root/doc/examples
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-06-30 10:11:55 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-06-30 10:11:55 -0400
commit4d4a23cd3e97bbc10d6307b828627c672097c5f5 (patch)
treeedccff8f095690ba489a55d4fee79bf6ba78b00d /doc/examples
parentb1741c1dc6c82737b69b20edfd817ef7e11ef28d (diff)
nearly complete page 6 / linear algebra + examples
fix the previous/next links
Diffstat (limited to 'doc/examples')
-rw-r--r--doc/examples/TutorialLinAlgComputeTwice.cpp23
-rw-r--r--doc/examples/TutorialLinAlgExComputeSolveError.cpp14
-rw-r--r--doc/examples/TutorialLinAlgExSolveColPivHouseholderQR.cpp6
-rw-r--r--doc/examples/TutorialLinAlgExSolveLDLT.cpp16
-rw-r--r--doc/examples/TutorialLinAlgInverseDeterminant.cpp16
-rw-r--r--doc/examples/TutorialLinAlgRankRevealing.cpp20
-rw-r--r--doc/examples/TutorialLinAlgSelfAdjointEigenSolver.cpp17
-rw-r--r--doc/examples/TutorialLinAlgSetThreshold.cpp19
8 files changed, 128 insertions, 3 deletions
diff --git a/doc/examples/TutorialLinAlgComputeTwice.cpp b/doc/examples/TutorialLinAlgComputeTwice.cpp
new file mode 100644
index 000000000..06ba6461a
--- /dev/null
+++ b/doc/examples/TutorialLinAlgComputeTwice.cpp
@@ -0,0 +1,23 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix2f A, b;
+ LLT<Matrix2f> llt;
+ A << 2, -1, -1, 3;
+ b << 1, 2, 3, 1;
+ cout << "Here is the matrix A:\n" << A << endl;
+ cout << "Here is the right hand side b:\n" << b << endl;
+ cout << "Computing LLT decomposition..." << endl;
+ llt.compute(A);
+ cout << "The solution is:\n" << llt.solve(b) << endl;
+ A(1,1)++;
+ cout << "The matrix A is now:\n" << A << endl;
+ cout << "Computing LLT decomposition..." << endl;
+ llt.compute(A);
+ cout << "The solution is now:\n" << llt.solve(b) << endl;
+}
diff --git a/doc/examples/TutorialLinAlgExComputeSolveError.cpp b/doc/examples/TutorialLinAlgExComputeSolveError.cpp
new file mode 100644
index 000000000..f362fb71a
--- /dev/null
+++ b/doc/examples/TutorialLinAlgExComputeSolveError.cpp
@@ -0,0 +1,14 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ MatrixXd A = MatrixXd::Random(100,100);
+ MatrixXd b = MatrixXd::Random(100,50);
+ MatrixXd x = A.fullPivLu().solve(b);
+ double relative_error = (A*x - b).norm() / b.norm(); // norm() is L2 norm
+ cout << "The relative error is:\n" << relative_error << endl;
+}
diff --git a/doc/examples/TutorialLinAlgExSolveColPivHouseholderQR.cpp b/doc/examples/TutorialLinAlgExSolveColPivHouseholderQR.cpp
index 29c22be41..3a99a94d7 100644
--- a/doc/examples/TutorialLinAlgExSolveColPivHouseholderQR.cpp
+++ b/doc/examples/TutorialLinAlgExSolveColPivHouseholderQR.cpp
@@ -10,8 +10,8 @@ int main()
Vector3f b;
A << 1,2,3, 4,5,6, 7,8,10;
b << 3, 3, 4;
- cout << "Here is the matrix A:" << endl << A << endl;
- cout << "Here is the vector b:" << endl << b << endl;
+ cout << "Here is the matrix A:\n" << A << endl;
+ cout << "Here is the vector b:\n" << b << endl;
Vector3f x = A.colPivHouseholderQr().solve(b);
- cout << "The solution is:" << endl << x << endl;
+ cout << "The solution is:\n" << x << endl;
}
diff --git a/doc/examples/TutorialLinAlgExSolveLDLT.cpp b/doc/examples/TutorialLinAlgExSolveLDLT.cpp
new file mode 100644
index 000000000..f8beacd27
--- /dev/null
+++ b/doc/examples/TutorialLinAlgExSolveLDLT.cpp
@@ -0,0 +1,16 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix2f A, b;
+ A << 2, -1, -1, 3;
+ b << 1, 2, 3, 1;
+ cout << "Here is the matrix A:\n" << A << endl;
+ cout << "Here is the right hand side b:\n" << b << endl;
+ Matrix2f x = A.ldlt().solve(b);
+ cout << "The solution is:\n" << x << endl;
+}
diff --git a/doc/examples/TutorialLinAlgInverseDeterminant.cpp b/doc/examples/TutorialLinAlgInverseDeterminant.cpp
new file mode 100644
index 000000000..43970ff05
--- /dev/null
+++ b/doc/examples/TutorialLinAlgInverseDeterminant.cpp
@@ -0,0 +1,16 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix3f A;
+ A << 1, 2, 1,
+ 2, 1, 0,
+ -1, 1, 2;
+ cout << "Here is the matrix A:\n" << A << endl;
+ cout << "The determinant of A is " << A.determinant() << endl;
+ cout << "The inverse of A is:\n" << A.inverse() << endl;
+} \ No newline at end of file
diff --git a/doc/examples/TutorialLinAlgRankRevealing.cpp b/doc/examples/TutorialLinAlgRankRevealing.cpp
new file mode 100644
index 000000000..c5165077f
--- /dev/null
+++ b/doc/examples/TutorialLinAlgRankRevealing.cpp
@@ -0,0 +1,20 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix3f A;
+ A << 1, 2, 5,
+ 2, 1, 4,
+ 3, 0, 3;
+ cout << "Here is the matrix A:\n" << A << endl;
+ FullPivLU<Matrix3f> lu_decomp(A);
+ cout << "The rank of A is " << lu_decomp.rank() << endl;
+ cout << "Here is a matrix whose columns form a basis of the null-space of A:\n"
+ << lu_decomp.kernel() << endl;
+ cout << "Here is a matrix whose columns form a basis of the column-space of A:\n"
+ << lu_decomp.image(A) << endl; // yes, have to pass the original A
+}
diff --git a/doc/examples/TutorialLinAlgSelfAdjointEigenSolver.cpp b/doc/examples/TutorialLinAlgSelfAdjointEigenSolver.cpp
new file mode 100644
index 000000000..e98444347
--- /dev/null
+++ b/doc/examples/TutorialLinAlgSelfAdjointEigenSolver.cpp
@@ -0,0 +1,17 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix2f A;
+ A << 1, 2, 2, 3;
+ cout << "Here is the matrix A:\n" << A << endl;
+ SelfAdjointEigenSolver<Matrix2f> eigensolver(A);
+ cout << "The eigenvalues of A are:\n" << eigensolver.eigenvalues() << endl;
+ cout << "Here's a matrix whose columns are eigenvectors of A "
+ << "corresponding to these eigenvalues:\n"
+ << eigensolver.eigenvectors() << endl;
+}
diff --git a/doc/examples/TutorialLinAlgSetThreshold.cpp b/doc/examples/TutorialLinAlgSetThreshold.cpp
new file mode 100644
index 000000000..e0927cf27
--- /dev/null
+++ b/doc/examples/TutorialLinAlgSetThreshold.cpp
@@ -0,0 +1,19 @@
+#include <iostream>
+#include <Eigen/Dense>
+
+using namespace std;
+using namespace Eigen;
+
+int main()
+{
+ Matrix2d A;
+ FullPivLU<Matrix2d> lu;
+ A << 2, 1,
+ 2, 0.9999999999;
+ lu.compute(A);
+ cout << "By default, the rank of A is found to be " << lu.rank() << endl;
+ cout << "Now recomputing the LU decomposition with threshold 1e-5" << endl;
+ lu.setThreshold(1e-5);
+ lu.compute(A);
+ cout << "The rank of A is found to be " << lu.rank() << endl;
+}