/home/coin/SVN-release/OS-2.4.0/Couenne/src/heuristics/CouenneFPcreateMILP.cpp

Go to the documentation of this file.
00001 /* $Id: CouenneFPcreateMILP.cpp 733 2011-07-03 22:32:52Z pbelotti $
00002  *
00003  * Name:    CouenneFPcreateMILP.cpp
00004  * Authors: Pietro Belotti
00005  *          Timo Berthold, ZIB Berlin
00006  * Purpose: create the MILP within the Feasibility Pump 
00007  * 
00008  * This file is licensed under the Eclipse Public License (EPL)
00009  */
00010 
00011 #include "CbcModel.hpp"
00012 
00013 #include "IpLapack.hpp"
00014 
00015 #include "CouenneSparseMatrix.hpp"
00016 #include "CouenneTNLP.hpp"
00017 #include "CouenneFeasPump.hpp"
00018 #include "CouenneProblem.hpp"
00019 #include "CouenneProblemElem.hpp"
00020 
00021 #define COUENNE_EIG_RATIO .1 // how much smaller than the largest eigenvalue should the minimum be set at?
00022 
00023 using namespace Couenne;
00024 
00026 void ComputeSquareRoot (const CouenneFeasPump *fp, CouenneSparseMatrix *hessian, CoinPackedVector *P);
00027 
00029 OsiSolverInterface *createCloneMILP (const CouenneFeasPump *fp, CbcModel *model, bool isMILP) {
00030 
00031   OsiSolverInterface *lp = model -> solver () -> clone ();
00032 
00033   // no data is available so far, retrieve it from the MILP solver
00034   // used as the linearization
00035 
00036   // Add q variables, each with coefficient 1 in the objective
00037 
00038   CoinPackedVector vec;
00039 
00040   for (int i = fp -> Problem () -> nVars (), j = 0; i--; ++j) {
00041 
00042     // column has to be added if:
00043     //
00044     // creating MIP AND (integer variable OR FP_DIST_ALL)
00045     // creating LP  AND fractional variable
00046 
00047     bool intVar = lp -> isInteger (j);
00048 
00049     if ((isMILP && (intVar || (fp -> compDistInt () == CouenneFeasPump::FP_DIST_ALL)))
00050         ||
00051         (!isMILP && !intVar))
00052       // (empty) coeff col vector, lb = 0, ub = inf, obj coeff
00053       lp -> addCol (vec, 0., COIN_DBL_MAX, 1.); 
00054   }
00055 
00056   // Set to zero all other variables' obj coefficient. This means we
00057   // just do it for the single variable in the reformulated
00058   // problem's linear relaxation (all other variables do not appear
00059   // in the objective)
00060 
00061   int objInd = fp -> Problem () -> Obj (0) -> Body () -> Index ();
00062 
00063   if (objInd >= 0)
00064     lp -> setObjCoeff (objInd, 0.);
00065 
00066   return lp;
00067 }
00068 
00069 
00071 void addDistanceConstraints (const CouenneFeasPump *fp, OsiSolverInterface *lp, double *sol, bool isMILP) {
00072 
00073   // Construct an (empty) Hessian. It will be modified later, but
00074   // the changes should be relatively easy for the case when
00075   // multHessMILP_ > 0 and there are no changes if multHessMILP_ == 0
00076 
00077   int n = fp -> Problem () -> nVars ();
00078 
00079   CoinPackedVector *P = new CoinPackedVector [n];
00080 
00081   // The MILP has to be changed the first time it is used.
00082   //
00083   // Suppose Ax >= b has m inequalities. In order to solve the
00084   // problem above, we need q new variables z_i and 2q inequalities
00085   //
00086   //   z_i >=   P^i (x - x^0)  or  P^i x - z_i <= P^i x^0 (*)
00087   //   z_i >= - P^i (x - x^0)                             (**)
00088   // 
00089   // (the latter being equivalent to
00090   //
00091   // - z_i <=   P^i (x - x^0)  or  P^i x + z_i >= P^i x^0 (***)
00092   //
00093   // so we'll use this instead as most coefficients don't change)
00094   // for each i, where q is the number of variables involved (either
00095   // q=|N|, the number of integer variables, or q=n, the number of
00096   // variables).
00097   //
00098   // We need to memorize the number of initial inequalities and of
00099   // variables, so that we know what (coefficients and rhs) to
00100   // change at every iteration.
00101 
00102   CoinPackedVector x0 (n, sol);
00103 
00104   // set objective coefficient if we are using a little Objective FP
00105 
00106   if (isMILP && (fp -> multObjFMILP () > 0.)) {
00107 
00108     int objInd = fp -> Problem () -> Obj (0) -> Body () -> Index ();
00109 
00110     if (objInd >= 0)
00111       lp -> setObjCoeff (objInd, fp -> multObjFMILP ());
00112   }
00113 
00114   if (isMILP && 
00115       (fp -> multHessMILP () > 0.) &&
00116       (fp -> nlp () -> optHessian ())) {
00117 
00118     // P is a convex combination, with weights multDistMILP_ and
00119     // multHessMILP_, of the distance and the Hessian respectively
00120 
00121     // obtain optHessian and compute its square root
00122 
00123     CouenneSparseMatrix *hessian = fp -> nlp () -> optHessian ();
00124 
00125     ComputeSquareRoot (fp, hessian, P);
00126 
00127   } else {
00128 
00129     // simply set P = I
00130 
00131     for (int i=0; i<n; i++)
00132       P[i].insert (i, 1.); 
00133   }
00134 
00135   // Add 2q inequalities
00136 
00137   for (int i = 0, j = n, k = j; k--; ++i) {
00138 
00139     // two rows have to be added if:
00140     //
00141     // amending MIP AND (integer variable OR FP_DIST_ALL)
00142     // amending ssLP  AND fractional variable
00143 
00144     bool intVar = lp -> isInteger (i);
00145 
00146     if ((isMILP && (intVar || (fp -> compDistInt () == CouenneFeasPump::FP_DIST_ALL)))
00147         ||
00148         (!isMILP && !intVar)) {
00149 
00150       // create vector with single entry of 1 at i-th position 
00151       CoinPackedVector &vec = P [i];
00152 
00153       if (vec.getNumElements () == 0)
00154         continue;
00155 
00156       // right-hand side equals <P^i,x^0>
00157       double PiX0 = sparseDotProduct (vec, x0); 
00158 
00159       // j is the index of the j-th extra variable z_j, used for z_j >=  P (x - x0)  ===> z_j - Px >= - Px_0 ==> -z_j + Px <= Px_0
00160       // Second inequality is                                    z_j >= -P (x - x0)                          ==>  z_j + Px >= Px_0
00161       vec.insert     (j,                         -1.); lp -> addRow (vec, -COIN_DBL_MAX,         PiX0); // (*)
00162       vec.setElement (vec.getNumElements () - 1, +1.); lp -> addRow (vec,          PiX0, COIN_DBL_MAX); // (***)
00163 
00164       ++j; // index of variable within problem (plus nVars_)
00165 
00166     } else if (intVar) { // implies (!isMILP)
00167 
00168       // fix integer variable to its value in iSol      
00169 
00170 #define INT_LP_BRACKET 0
00171 
00172       lp -> setColLower (i, sol [i] - INT_LP_BRACKET);
00173       lp -> setColUpper (i, sol [i] + INT_LP_BRACKET);
00174     }
00175   }
00176 
00177   delete [] P;
00178 }
00179 
00180 
00181 #define GRADIENT_WEIGHT 1
00182 
00183 void ComputeSquareRoot (const CouenneFeasPump *fp, 
00184                         CouenneSparseMatrix *hessian, 
00185                         CoinPackedVector *P) {
00186   int 
00187     objInd = fp -> Problem () -> Obj (0) -> Body () -> Index (),
00188     n      = fp -> Problem () -> nVars ();
00189 
00190   //assert (objInd >= 0);
00191 
00192   double *val = hessian -> val ();
00193   int    *row = hessian -> row ();
00194   int    *col = hessian -> col ();
00195   int     num = hessian -> num ();
00196 
00197   //printf ("compute square root:\n");
00198 
00199   // Remove objective's row and column (equivalent to taking the
00200   // Lagrangian's Hessian, setting f(x) = x_z = c, and recomputing the
00201   // hessian). 
00202 
00203   double maxElem = 0.; // used in adding diagonal element of x_z
00204 
00205   for (int i=num; i--; ++row, ++col, ++val) {
00206 
00207     //printf ("elem: %d, %d --> %g\n", *row, *col, *val);
00208 
00209     if ((*row == objInd) || 
00210         (*col == objInd))
00211 
00212       *val = 0;
00213 
00214     else if (fabs (*val) > maxElem)
00215       maxElem = fabs (*val);
00216   }
00217 
00218   val -= num;
00219   row -= num;
00220   col -= num;
00221 
00222   // fill an input to Lapack/Blas procedures using hessian
00223 
00224   double *A = (double *) malloc (n*n * sizeof (double));
00225 
00226   CoinZeroN (A, n*n);
00227 
00228   // Add Hessian part -- only lower triangular part
00229   for (int i=0; i<num; ++i, ++row, ++col, ++val)
00230     if (*col <= *row)
00231       A [*col * n + *row] = fp -> multHessMILP () * *val;
00232 
00233   // Add distance part
00234   for (int i=0; i<n; ++i)
00235     A [i * (n+1)] += fp -> multDistMILP ();
00236 
00237   // Add gradient-parallel term to the hessian, (x_z - x_z_0)^2. This
00238   // amounts to setting the diagonal element to GRADIENT_WEIGHT. Don't
00239   // do it directly on hessian
00240 
00241   if (objInd >= 0)
00242     A [objInd * (n+1)] = maxElem * GRADIENT_WEIGHT * n;
00243 
00244   // call Lapack/Blas routines
00245   double *eigenval = (double *) malloc (n   * sizeof (double));
00246   int status;
00247 
00248   // compute eigenvalues and eigenvectors
00249   Ipopt::IpLapackDsyev (true, n, A, n, eigenval, status);
00250 
00251   if      (status < 0) printf ("Couenne: warning, argument %d illegal\n",                     -status);
00252   else if (status > 0) printf ("Couenne: warning, dsyev did not converge (error code: %d)\n",  status);
00253 
00254   // define a new matrix B = E' * D, where E' is E transposed,
00255   //
00256   // E = eigenvector matrix
00257   // D = diagonal with square roots of eigenvalues
00258   //
00259   // Eventually, the square root is given by E' D E
00260 
00261   double *B = (double *) malloc (n*n * sizeof(double));
00262 
00263   double *eigenvec = A; // as overwritten by dsyev;
00264 
00265   // 
00266   // eigenvec is column major, hence post-multiplying it by D equals
00267   // multiplying each column i by the i-th eigenvalue
00268   //
00269 
00270   // if all eigenvalues are nonpositive, set them all to one
00271 
00272   double
00273     MinEigVal = eigenval [0],
00274     MaxEigVal = eigenval [n-1];
00275 
00276   for (int j=1; j<n; j++)
00277     assert (eigenval [j-1] <= eigenval [j]);
00278 
00279   if (MaxEigVal <= 0.)
00280 
00281     // in this case it makes sense to invert each eigenvalue
00282     // (i.e. take its inverse) and change its sign, as the steepest
00283     // descent should correspond to the thinnest direction
00284 
00285     for (int j=0; j<n; j++)
00286       eigenval [j] = 1. / (.1 - eigenval [j]);
00287 
00288   else {
00289 
00290     // set all not-too-positive ones to a fraction of the maximum
00291     // ("un-thins" the level curves defined by the HL)
00292 
00293     MinEigVal = MaxEigVal * COUENNE_EIG_RATIO;
00294 
00295     if (eigenval [0] <= MinEigVal) 
00296       for (int j=0; eigenval [j] <= MinEigVal; j++)
00297         eigenval [j] = MinEigVal;
00298   }
00299 
00300   // Now obtain sqrt (A)
00301 
00302   for (int j=0; j<n; ++j) {
00303 
00304     register double sqrtEig = sqrt (eigenval [j]);
00305 
00306     for (int i=n; i--;)
00307       *B++ = sqrtEig * eigenvec [i*n+j];
00308   }
00309 
00310   B -= n*n;
00311 
00312   // TODO: set up B as    row-major sparse matrix and 
00313   //              E as column-major sparse matrix
00314   //
00315   // Otherwise this multiplication is O(n^3)
00316 
00317   // Now compute B * E
00318 
00319   for     (int i=0; i<n; ++i)
00320     for   (int j=0; j<n; ++j) { 
00321 
00322       // multiply i-th row of B by j-th column of E
00323 
00324       double elem = 0.;
00325 
00326       for (int k=0; k<n; ++k)
00327         elem += B [i + k * n] * eigenvec [j * n + k];
00328 
00329       if (fabs (elem) > COUENNE_EPS)
00330         P [i]. insert (j, elem);
00331     }
00332 
00333   if (fp -> Problem () -> Jnlst () -> ProduceOutput (Ipopt::J_STRONGWARNING, J_NLPHEURISTIC)) {
00334 
00335     printf ("P:\n");
00336 
00337 
00338     printf ("P^{1/2}:\n");
00339 
00340   }
00341 
00342   free (eigenval);
00343   free (A);
00344   free (B);
00345 }

Generated on Thu Sep 22 03:05:58 2011 by  doxygen 1.4.7