source: pacpussensors/trunk/Vislab/lib3dv-1.2.0/lib3dv/eigen/doc/examples/matrixfree_cg.cpp

Last change on this file was 136, checked in by ldecherf, 8 years ago

Doc

File size: 5.4 KB
Line 
1#include <iostream>
2#include <Eigen/Core>
3#include <Eigen/Dense>
4#include <Eigen/IterativeLinearSolvers>
5
6class MatrixReplacement;
7template<typename Rhs> class MatrixReplacement_ProductReturnType;
8
9namespace Eigen {
10namespace internal {
11 template<>
12 struct traits<MatrixReplacement> : Eigen::internal::traits<Eigen::SparseMatrix<double> >
13 {};
14
15 template <typename Rhs>
16 struct traits<MatrixReplacement_ProductReturnType<Rhs> > {
17 // The equivalent plain objet type of the product. This type is used if the product needs to be evaluated into a temporary.
18 typedef Eigen::Matrix<typename Rhs::Scalar, Eigen::Dynamic, Rhs::ColsAtCompileTime> ReturnType;
19 };
20}
21}
22
23// Inheriting EigenBase should not be needed in the future.
24class MatrixReplacement : public Eigen::EigenBase<MatrixReplacement> {
25public:
26 // Expose some compile-time information to Eigen:
27 typedef double Scalar;
28 typedef double RealScalar;
29 enum {
30 ColsAtCompileTime = Eigen::Dynamic,
31 RowsAtCompileTime = Eigen::Dynamic,
32 MaxColsAtCompileTime = Eigen::Dynamic,
33 MaxRowsAtCompileTime = Eigen::Dynamic
34 };
35
36 Index rows() const { return 4; }
37 Index cols() const { return 4; }
38
39 void resize(Index a_rows, Index a_cols)
40 {
41 // This method should not be needed in the future.
42 assert(a_rows==0 && a_cols==0 || a_rows==rows() && a_cols==cols());
43 }
44
45 // In the future, the return type should be Eigen::Product<MatrixReplacement,Rhs>
46 template<typename Rhs>
47 MatrixReplacement_ProductReturnType<Rhs> operator*(const Eigen::MatrixBase<Rhs>& x) const {
48 return MatrixReplacement_ProductReturnType<Rhs>(*this, x.derived());
49 }
50
51};
52
53// The proxy class representing the product of a MatrixReplacement with a MatrixBase<>
54template<typename Rhs>
55class MatrixReplacement_ProductReturnType : public Eigen::ReturnByValue<MatrixReplacement_ProductReturnType<Rhs> > {
56public:
57 typedef MatrixReplacement::Index Index;
58
59 // The ctor store references to the matrix and right-hand-side object (usually a vector).
60 MatrixReplacement_ProductReturnType(const MatrixReplacement& matrix, const Rhs& rhs)
61 : m_matrix(matrix), m_rhs(rhs)
62 {}
63
64 Index rows() const { return m_matrix.rows(); }
65 Index cols() const { return m_rhs.cols(); }
66
67 // This function is automatically called by Eigen. It must evaluate the product of matrix * rhs into y.
68 template<typename Dest>
69 void evalTo(Dest& y) const
70 {
71 y.setZero(4);
72
73 y(0) += 2 * m_rhs(0); y(1) += 1 * m_rhs(0);
74 y(0) += 1 * m_rhs(1); y(1) += 2 * m_rhs(1); y(2) += 1 * m_rhs(1);
75 y(1) += 1 * m_rhs(2); y(2) += 2 * m_rhs(2); y(3) += 1 * m_rhs(2);
76 y(2) += 1 * m_rhs(3); y(3) += 2 * m_rhs(3);
77 }
78
79protected:
80 const MatrixReplacement& m_matrix;
81 typename Rhs::Nested m_rhs;
82};
83
84
85/*****/
86
87// This class simply warp a diagonal matrix as a Jacobi preconditioner.
88// In the future such simple and generic wrapper should be shipped within Eigen itsel.
89template <typename _Scalar>
90class MyJacobiPreconditioner
91{
92 typedef _Scalar Scalar;
93 typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> Vector;
94 typedef typename Vector::Index Index;
95
96 public:
97 // this typedef is only to export the scalar type and compile-time dimensions to solve_retval
98 typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixType;
99
100 MyJacobiPreconditioner() : m_isInitialized(false) {}
101
102 void setInvDiag(const Eigen::VectorXd &invdiag) {
103 m_invdiag=invdiag;
104 m_isInitialized=true;
105 }
106
107 Index rows() const { return m_invdiag.size(); }
108 Index cols() const { return m_invdiag.size(); }
109
110 template<typename MatType>
111 MyJacobiPreconditioner& analyzePattern(const MatType& ) { return *this; }
112
113 template<typename MatType>
114 MyJacobiPreconditioner& factorize(const MatType& mat) { return *this; }
115
116 template<typename MatType>
117 MyJacobiPreconditioner& compute(const MatType& mat) { return *this; }
118
119 template<typename Rhs, typename Dest>
120 void _solve(const Rhs& b, Dest& x) const
121 {
122 x = m_invdiag.array() * b.array() ;
123 }
124
125 template<typename Rhs> inline const Eigen::internal::solve_retval<MyJacobiPreconditioner, Rhs>
126 solve(const Eigen::MatrixBase<Rhs>& b) const
127 {
128 eigen_assert(m_isInitialized && "MyJacobiPreconditioner is not initialized.");
129 eigen_assert(m_invdiag.size()==b.rows()
130 && "MyJacobiPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
131 return Eigen::internal::solve_retval<MyJacobiPreconditioner, Rhs>(*this, b.derived());
132 }
133
134 protected:
135 Vector m_invdiag;
136 bool m_isInitialized;
137};
138
139namespace Eigen {
140namespace internal {
141
142template<typename _MatrixType, typename Rhs>
143struct solve_retval<MyJacobiPreconditioner<_MatrixType>, Rhs>
144 : solve_retval_base<MyJacobiPreconditioner<_MatrixType>, Rhs>
145{
146 typedef MyJacobiPreconditioner<_MatrixType> Dec;
147 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
148
149 template<typename Dest> void evalTo(Dest& dst) const
150 {
151 dec()._solve(rhs(),dst);
152 }
153};
154
155}
156}
157
158
159/*****/
160
161
162int main()
163{
164 MatrixReplacement A;
165 Eigen::VectorXd b(4), x;
166 b << 1, 1, 1, 1;
167
168 // solve Ax = b using CG with matrix-free version:
169 Eigen::ConjugateGradient < MatrixReplacement, Eigen::Lower|Eigen::Upper, MyJacobiPreconditioner<double> > cg;
170
171 Eigen::VectorXd invdiag(4);
172 invdiag << 1./3., 1./4., 1./4., 1./3.;
173
174 cg.preconditioner().setInvDiag(invdiag);
175 cg.compute(A);
176 x = cg.solve(b);
177
178 std::cout << "#iterations: " << cg.iterations() << std::endl;
179 std::cout << "estimated error: " << cg.error() << std::endl;
180}
Note: See TracBrowser for help on using the repository browser.