ConstrainedConjGrad.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10/* NOTE The functions of this file have been adapted from the GMM++ library */
11
12//========================================================================
13//
14// Copyright (C) 2002-2007 Yves Renard
15//
16// This file is a part of GETFEM++
17//
18// Getfem++ is free software; you can redistribute it and/or modify
19// it under the terms of the GNU Lesser General Public License as
20// published by the Free Software Foundation; version 2.1 of the License.
21//
22// This program is distributed in the hope that it will be useful,
23// but WITHOUT ANY WARRANTY; without even the implied warranty of
24// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25// GNU Lesser General Public License for more details.
26// You should have received a copy of the GNU Lesser General Public
27// License along with this program; if not, write to the Free Software
28// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301,
29// USA.
30//
31//========================================================================
32
33#include "../../../../Eigen/src/Core/util/NonMPL2.h"
34
35#ifndef EIGEN_CONSTRAINEDCG_H
36#define EIGEN_CONSTRAINEDCG_H
37
38#include <Eigen/Core>
39
40namespace Eigen {
41
42namespace internal {
43
50template <typename CMatrix, typename CINVMatrix>
51void pseudo_inverse(const CMatrix &C, CINVMatrix &CINV)
52{
53 // optimisable : copie de la ligne, precalcul de C * trans(C).
54 typedef typename CMatrix::Scalar Scalar;
55 typedef typename CMatrix::Index Index;
56 // FIXME use sparse vectors ?
57 typedef Matrix<Scalar,Dynamic,1> TmpVec;
58
59 Index rows = C.rows(), cols = C.cols();
60
61 TmpVec d(rows), e(rows), l(cols), p(rows), q(rows), r(rows);
62 Scalar rho, rho_1, alpha;
63 d.setZero();
64
65 typedef Triplet<double> T;
66 std::vector<T> tripletList;
67
68 for (Index i = 0; i < rows; ++i)
69 {
70 d[i] = 1.0;
71 rho = 1.0;
72 e.setZero();
73 r = d;
74 p = d;
75
76 while (rho >= 1e-38)
77 { /* conjugate gradient to compute e */
78 /* which is the i-th row of inv(C * trans(C)) */
79 l = C.transpose() * p;
80 q = C * l;
81 alpha = rho / p.dot(q);
82 e += alpha * p;
83 r += -alpha * q;
84 rho_1 = rho;
85 rho = r.dot(r);
86 p = (rho/rho_1) * p + r;
87 }
88
89 l = C.transpose() * e; // l is the i-th row of CINV
90 // FIXME add a generic "prune/filter" expression for both dense and sparse object to sparse
91 for (Index j=0; j<l.size(); ++j)
92 if (l[j]<1e-15)
93 tripletList.push_back(T(i,j,l(j)));
94
95
96 d[i] = 0.0;
97 }
98 CINV.setFromTriplets(tripletList.begin(), tripletList.end());
99}
100
101
102
108template<typename TMatrix, typename CMatrix,
109 typename VectorX, typename VectorB, typename VectorF>
110void constrained_cg(const TMatrix& A, const CMatrix& C, VectorX& x,
111 const VectorB& b, const VectorF& f, IterationController &iter)
112{
113 using std::sqrt;
114 typedef typename TMatrix::Scalar Scalar;
115 typedef typename TMatrix::Index Index;
116 typedef Matrix<Scalar,Dynamic,1> TmpVec;
117
118 Scalar rho = 1.0, rho_1, lambda, gamma;
119 Index xSize = x.size();
120 TmpVec p(xSize), q(xSize), q2(xSize),
121 r(xSize), old_z(xSize), z(xSize),
122 memox(xSize);
123 std::vector<bool> satured(C.rows());
124 p.setZero();
125 iter.setRhsNorm(sqrt(b.dot(b))); // gael vect_sp(PS, b, b)
126 if (iter.rhsNorm() == 0.0) iter.setRhsNorm(1.0);
127
128 SparseMatrix<Scalar,RowMajor> CINV(C.rows(), C.cols());
129 pseudo_inverse(C, CINV);
130
131 while(true)
132 {
133 // computation of residual
134 old_z = z;
135 memox = x;
136 r = b;
137 r += A * -x;
138 z = r;
139 bool transition = false;
140 for (Index i = 0; i < C.rows(); ++i)
141 {
142 Scalar al = C.row(i).dot(x) - f.coeff(i);
143 if (al >= -1.0E-15)
144 {
145 if (!satured[i])
146 {
147 satured[i] = true;
148 transition = true;
149 }
150 Scalar bb = CINV.row(i).dot(z);
151 if (bb > 0.0)
152 // FIXME: we should allow that: z += -bb * C.row(i);
153 for (typename CMatrix::InnerIterator it(C,i); it; ++it)
154 z.coeffRef(it.index()) -= bb*it.value();
155 }
156 else
157 satured[i] = false;
158 }
159
160 // descent direction
161 rho_1 = rho;
162 rho = r.dot(z);
163
164 if (iter.finished(rho)) break;
165
166 if (iter.noiseLevel() > 0 && transition) std::cerr << "CCG: transition\n";
167 if (transition || iter.first()) gamma = 0.0;
168 else gamma = (std::max)(0.0, (rho - old_z.dot(z)) / rho_1);
169 p = z + gamma*p;
170
171 ++iter;
172 // one dimensionnal optimization
173 q = A * p;
174 lambda = rho / q.dot(p);
175 for (Index i = 0; i < C.rows(); ++i)
176 {
177 if (!satured[i])
178 {
179 Scalar bb = C.row(i).dot(p) - f[i];
180 if (bb > 0.0)
181 lambda = (std::min)(lambda, (f.coeff(i)-C.row(i).dot(x)) / bb);
182 }
183 }
184 x += lambda * p;
185 memox -= x;
186 }
187}
188
189} // end namespace internal
190
191} // end namespace Eigen
192
193#endif // EIGEN_CONSTRAINEDCG_H