ConstrainedConjGrad.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 /* NOTE The functions of this file have been adapted from the GMM++ library */
00026 
00027 //========================================================================
00028 //
00029 // Copyright (C) 2002-2007 Yves Renard
00030 //
00031 // This file is a part of GETFEM++
00032 //
00033 // Getfem++ is free software; you can redistribute it and/or modify
00034 // it under the terms of the GNU Lesser General Public License as
00035 // published by the Free Software Foundation; version 2.1 of the License.
00036 //
00037 // This program is distributed in the hope that it will be useful,
00038 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00039 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00040 // GNU Lesser General Public License for more details.
00041 // You should have received a copy of the GNU Lesser General Public
00042 // License along with this program; if not, write to the Free Software
00043 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301,
00044 // USA.
00045 //
00046 //========================================================================
00047 
00048 #ifndef EIGEN_CONSTRAINEDCG_H
00049 #define EIGEN_CONSTRAINEDCG_H
00050 
00051 #include <Eigen/Core>
00052 
00053 namespace internal {
00054 
00061 template <typename CMatrix, typename CINVMatrix>
00062 void pseudo_inverse(const CMatrix &C, CINVMatrix &CINV)
00063 {
00064   // optimisable : copie de la ligne, precalcul de C * trans(C).
00065   typedef typename CMatrix::Scalar Scalar;
00066   typedef typename CMatrix::Index Index;
00067   // FIXME use sparse vectors ?
00068   typedef Matrix<Scalar,Dynamic,1> TmpVec;
00069 
00070   Index rows = C.rows(), cols = C.cols();
00071 
00072   TmpVec d(rows), e(rows), l(cols), p(rows), q(rows), r(rows);
00073   Scalar rho, rho_1, alpha;
00074   d.setZero();
00075 
00076   CINV.startFill(); // FIXME estimate the number of non-zeros
00077   for (Index i = 0; i < rows; ++i)
00078   {
00079     d[i] = 1.0;
00080     rho = 1.0;
00081     e.setZero();
00082     r = d;
00083     p = d;
00084 
00085     while (rho >= 1e-38)
00086     { /* conjugate gradient to compute e             */
00087       /* which is the i-th row of inv(C * trans(C))  */
00088       l = C.transpose() * p;
00089       q = C * l;
00090       alpha = rho / p.dot(q);
00091       e +=  alpha * p;
00092       r += -alpha * q;
00093       rho_1 = rho;
00094       rho = r.dot(r);
00095       p = (rho/rho_1) * p + r;
00096     }
00097 
00098     l = C.transpose() * e; // l is the i-th row of CINV
00099     // FIXME add a generic "prune/filter" expression for both dense and sparse object to sparse
00100     for (Index j=0; j<l.size(); ++j)
00101       if (l[j]<1e-15)
00102         CINV.fill(i,j) = l[j];
00103 
00104     d[i] = 0.0;
00105   }
00106   CINV.endFill();
00107 }
00108 
00109 
00110 
00116 template<typename TMatrix, typename CMatrix,
00117          typename VectorX, typename VectorB, typename VectorF>
00118 void constrained_cg(const TMatrix& A, const CMatrix& C, VectorX& x,
00119                        const VectorB& b, const VectorF& f, IterationController &iter)
00120 {
00121   typedef typename TMatrix::Scalar Scalar;
00122   typedef typename TMatrix::Index Index;
00123   typedef Matrix<Scalar,Dynamic,1>  TmpVec;
00124 
00125   Scalar rho = 1.0, rho_1, lambda, gamma;
00126   Index xSize = x.size();
00127   TmpVec  p(xSize), q(xSize), q2(xSize),
00128           r(xSize), old_z(xSize), z(xSize),
00129           memox(xSize);
00130   std::vector<bool> satured(C.rows());
00131   p.setZero();
00132   iter.setRhsNorm(sqrt(b.dot(b))); // gael vect_sp(PS, b, b)
00133   if (iter.rhsNorm() == 0.0) iter.setRhsNorm(1.0);
00134 
00135   SparseMatrix<Scalar,RowMajor> CINV(C.rows(), C.cols());
00136   pseudo_inverse(C, CINV);
00137 
00138   while(true)
00139   {
00140     // computation of residual
00141     old_z = z;
00142     memox = x;
00143     r = b;
00144     r += A * -x;
00145     z = r;
00146     bool transition = false;
00147     for (Index i = 0; i < C.rows(); ++i)
00148     {
00149       Scalar al = C.row(i).dot(x) - f.coeff(i);
00150       if (al >= -1.0E-15)
00151       {
00152         if (!satured[i])
00153         {
00154           satured[i] = true;
00155           transition = true;
00156         }
00157         Scalar bb = CINV.row(i).dot(z);
00158         if (bb > 0.0)
00159           // FIXME: we should allow that: z += -bb * C.row(i);
00160           for (typename CMatrix::InnerIterator it(C,i); it; ++it)
00161             z.coeffRef(it.index()) -= bb*it.value();
00162       }
00163       else
00164         satured[i] = false;
00165     }
00166 
00167     // descent direction
00168     rho_1 = rho;
00169     rho = r.dot(z);
00170 
00171     if (iter.finished(rho)) break;
00172 
00173     if (iter.noiseLevel() > 0 && transition) std::cerr << "CCG: transition\n";
00174     if (transition || iter.first()) gamma = 0.0;
00175     else gamma = std::max(0.0, (rho - old_z.dot(z)) / rho_1);
00176     p = z + gamma*p;
00177 
00178     ++iter;
00179     // one dimensionnal optimization
00180     q = A * p;
00181     lambda = rho / q.dot(p);
00182     for (Index i = 0; i < C.rows(); ++i)
00183     {
00184       if (!satured[i])
00185       {
00186         Scalar bb = C.row(i).dot(p) - f[i];
00187         if (bb > 0.0)
00188           lambda = std::min(lambda, (f.coeff(i)-C.row(i).dot(x)) / bb);
00189       }
00190     }
00191     x += lambda * p;
00192     memox -= x;
00193   }
00194 }
00195 
00196 } // end namespace internal
00197 
00198 #endif // EIGEN_CONSTRAINEDCG_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:30:59