• Main Page
  • Related Pages
  • Modules
  • Data Structures
  • Files
  • File List
  • Globals

sst/elements/genericProc/programs/MTGL/mtgl/ufl.h

00001 /*  _________________________________________________________________________ 
00002  * 
00003  *  MTGL: The MultiThreaded Graph Library
00004  *  Copyright (c) 2008 Sandia Corporation. 
00005  *  This software is distributed under the BSD License. 
00006  *  Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 
00007  *  the U.S. Government retains certain rights in this software. 
00008  *  For more information, see the README file in the top MTGL directory. 
00009  *  _________________________________________________________________________ 
00010  */ 
00011 
00012 // This code is licensed under the Common Public License
00013 // It originated from COIN-OR
00014 
00015 // Copyright (C) 2000, International Business Machines
00016 // Corporation and others.  All Rights Reserved.
00017 
00018 #ifndef MTGL_UFL_H
00019 #define MTGL_UFL_H
00020 
00021 #define PARCUTOFF 5000
00022 
00023 #include <vector>       // including sys/times.h before this led to compiler
00024                         // errors on Linux
00025 #include <cstdio>
00026 #include <cfloat>
00027 #include <string>
00028 #include <map>
00029 #include <cstring>
00030 
00031 #include <mtgl/util.hpp>
00032 #include <mtgl/quicksort.hpp>
00033 #include <mtgl/rand_fill.hpp>
00034 #include <mtgl/SMVkernel.h>
00035 #include <mtgl/VolVolume.h>
00036 
00037 #ifndef WIN32
00038 #include <sys/times.h>
00039 #endif
00040 
00041 #define ALWAYS_SERVE 1
00042 #define DONT_ALWAYS_SERVE 0
00043 
00044 using namespace mtgl;
00045 
00046 // parameters controlled by the user
00047 class UFL_parms {
00048 public:
00049   std::string fdata; // file with the data
00050   std::string dualfile; // file with an initial dual solution
00051   std::string dual_savefile; // file to save final dual solution
00052   std::string int_savefile; // file to save primal integer solution
00053   int h_iter; // number of times that the primal heuristic will be run
00054   // after termination of the volume algorithm
00055    
00056   UFL_parms(): fdata(""), dualfile(""), dual_savefile(""), int_savefile("sensors.txt"), h_iter(100){}
00057 
00058   UFL_parms(const char* filename);
00059   ~UFL_parms() {}
00060 };
00061 
00062 
00063 
00064 template <typename vtype, typename adjlist_t>
00065 class UFL : public VOL_user_hooks {
00066 public:
00067   typedef typename adjlist_t::column_iterator column_iterator;
00068   typedef typename adjlist_t::value_iterator value_iterator;
00069   // for all hooks: return value of -1 means that volume should quit
00070   // compute reduced costs
00071    int compute_rc(const VOL_dvector& u, VOL_dvector& rc);
00072   // solve relaxed problem
00073    int solve_subproblem(const VOL_dvector& u, const VOL_dvector& rc,
00074                         double& lcost, VOL_dvector& x, VOL_dvector&v,
00075                         double& pcost);
00076    //DEPRECATED
00077 /*    int solve_subproblem_with_mask(const VOL_dvector& u, const VOL_dvector& rc, */
00078 /*                      double& lcost, VOL_dvector& x, VOL_dvector&v, */
00079 /*                                double& pcost, int **load_balance_masks = 0, */
00080 /*                                int *load_balance_totals = 0, */
00081 /*                                int **load_balance_parallelize_flags = 0,     */
00082 /*                                int num_load_balance_masks = 0); */
00083 
00084    void process_opening_serving_of_nodes(int mask_iteration, int mask_vertex_index, 
00085                                          int &num_open, double &sum_fcosts, double &value, 
00086                                          VOL_ivector &sol, const double *rdist,
00087                                          int **load_balance_masks, int *load_balance_totals, 
00088                                          int **load_balance_parallelize_flags, int num_load_balance_masks);
00089 
00090   // primal heuristic
00091   // return DBL_MAX in heur_val if feas sol wasn't/was found 
00092    int heuristics(const VOL_problem& p, 
00093                   const VOL_dvector& x, double& heur_val, double lb);
00094 
00095    void init_u(VOL_dvector& u);
00096 
00097    double open_more_servers(VOL_ivector& nsol,
00098                             int *unserved, int num_unserved);
00099    double assign_service(VOL_ivector& nsol, const VOL_dvector& x);
00100    int * find_unserved(VOL_ivector& nsol, int num_servers, int& num_unserved);
00101    void separate_by_degree();
00102    void separate_by_degree(int *to_sep, int num_to_sep,
00103                                          int *&hd, int& num_hd,
00104                                          int *&ld, int& num_ld);
00105 
00106   // original data for uncapacitated facility location
00107 public: 
00108   VOL_dvector fcost; // cost for opening facilities
00109   VOL_dvector dist; // cost for connecting a customer to a facility
00110   adjlist_t& _dist;
00111   std::vector<adjlist_t > _sidecon; // Additional side constr.
00112   std::vector<vtype> _ub; // Upper bounds for side constraints
00113   VOL_ivector fix; // vector saying if some variables should be fixed
00114   // if fix=-1 nothing is fixed
00115   int ncust, nloc; // number of customers, number of locations
00116   int nnz;      // number of nonzero (customer, facility) pairs
00117   int _p;          // _p>0 for p-median problems 
00118   VOL_ivector ix;   // best integer feasible solution so far
00119   double      icost;  // value of best integer feasible solution 
00120   double      gap;    // percentage (icost-lb)/lb
00121   double     *init_lagrange;  // user may want to initialize multipliers
00122   //int        *high_degree;    // indices of high degree locations
00123   //int        *low_degree;     // indices of low degree locations
00124   //int        n_hd;          // number of high degree locations
00125   //int        n_ld;          // number of low degree locations
00126   
00127   int        twenty_fifth;   //should take this out once we decide where it should go
00128 
00129 public:
00130   UFL(adjlist_t& d, double *init_l=0) : _dist(d), icost(DBL_MAX),
00131                                         init_lagrange(init_l)//,
00132                                         //n_hd(0),n_ld(0),
00133                                         //high_degree(0), low_degree(0), 
00134                                         //twenty_fifth(0)
00135   {
00136   }
00137   virtual ~UFL() 
00138   {
00139         //delete[] high_degree;
00140         //delete[] low_degree;
00141   }  
00142 
00143   //int get_n_hd() const { return n_hd; }
00144   //int get_n_ld() const { return n_ld; }
00145   //int* get_high_degree() const { return high_degree; }
00146   //int* get_low_degree()  const { return low_degree;  }
00147 };
00148 
00149 //#############################################################################
00150 //########  Member functions          #########################################
00151 //#############################################################################
00152 
00153 
00154 extern int newDummyID;
00155 extern map<int, int> newID;
00156 extern map<int, int>::iterator newIDiter;
00157 extern vector<int> origID;
00158 extern double gap;      // e.g. 0.1 --> stop when within 10% of optimal
00159 
00160 char *UFL_preprocess_data(const char* fname,  int num_facilities);
00161 
00162 template <typename vtype, typename adjlist>
00163 int ufl_module(int nloc, int ncust, int nnz,
00164            int& num_facilities_to_open, double *opening_costs, 
00165            adjlist& service_costs, double gap,
00166            bool* open_facilities, int *server,
00167            double *frac_solution=0,
00168            double *init_lagrange=0,
00169            int **load_balance_masks = 0,
00170            int *load_balance_totals = 0,
00171            int **load_balance_parallelize_flags = 0,    
00172            int num_load_balance_masks = 0
00173 )
00174 {
00175    typedef typename adjlist::column_iterator column_iterator;
00176    typedef typename adjlist::value_iterator  value_iterator;
00177 
00178 
00179    printf("sizeof(volp): %lu\n", (long unsigned)sizeof(VOL_problem));
00180    //UFL_parms ufl_par("ufl.par"); // for debugging
00181    UFL_parms ufl_par;   // SPOT - omit parameter file
00182    UFL<vtype, adjlist>  ufl_data(service_costs, init_lagrange);
00183    int num_sidecon = 0; // Default is no side constraints
00184 
00185    newDummyID = nloc;   // dummy concept irrelevant in this context
00186    // read data 
00187    printf("ufl_module: nloc: %d\n", nloc);
00188    ufl_data.nloc = nloc;
00189    ufl_data.ncust = ncust;
00190    ufl_data.nnz = nnz;
00191    ufl_data.gap = gap;
00192    VOL_dvector& fcost = ufl_data.fcost;
00193    fcost.allocate(nloc);
00194 
00195    if (opening_costs) {
00196         #pragma mta assert nodep
00197         for (int i=0; i<nloc; i++) {
00198                 fcost[i] = opening_costs[i];
00199         }
00200         ufl_data._p = 0;
00201    } else {
00202         ufl_data._p = num_facilities_to_open;
00203         #pragma mta assert nodep
00204         for (int i=0; i<nloc; i++) {
00205                 fcost[i] = 0.0;
00206         }
00207    }
00208    //ufl_data.separate_by_degree();
00209 
00210    // create the VOL_problem from the parameter file
00211    //VOL_problem volp("ufl.par"); // for debugging
00212    VOL_problem volp;    // SPOT - omit parameter file
00213    // volp.psize = ufl_data.nloc + ufl_data.nloc*ufl_data.ncust;
00214    volp.psize = ufl_data.nloc + ufl_data.nnz;
00215    volp.dsize = ufl_data.ncust + num_sidecon;
00216    bool ifdual = false;
00217    int state = 12345;
00218    if (ufl_par.dualfile.length() > 0) {
00219      // read dual solution
00220       ifdual = true;
00221       VOL_dvector& dinit = volp.dsol;
00222       dinit.allocate(volp.dsize);
00223 
00224       FILE * file = fopen(ufl_par.dualfile.c_str(), "r");
00225       if (!file) {
00226          printf("Failure to open file: %s\n ", ufl_par.dualfile.c_str());
00227          abort();
00228       }
00229       const int dsize = volp.dsize;
00230       int idummy;
00231       for (int i = 0; i < dsize; ++i) {
00232          int rc = fscanf(file, "%i%lf", &idummy, &dinit[i]);
00233          if (rc != 2) {
00234            printf("ERROR: ufl_module: failed to scan all fields\n");
00235          }
00236       }
00237       fclose(file);
00238    } else {
00239         // initialize lagrange multipliers in a slightly random way to 
00240         // avoid the cycling problem associated with self-loops
00241         VOL_dvector& dinit = volp.dsol;  
00242         dinit.allocate(volp.dsize);
00243         int state = 12345;
00244         int vdsize = volp.dsize;
00245         double r;
00246 
00247   random_value* rand =
00248     (random_value*) malloc(volp.dsize * sizeof(random_value));
00249         rand_fill::generate(volp.dsize, rand);
00250 
00251         #pragma mta assert nodep
00252         for (int i = 0; i < vdsize; ++i) {
00253                 r= rand[i] / (double) INT_MAX;
00254                 dinit[i] = r;
00255         }
00256    }
00257 
00258    // We set dual bounds for side constraints. Assume primal <= constraints.
00259 
00260    // This would be the right place to set bounds on the dual variables
00261    // For UFL all the relaxed constraints are equalities, so the bounds are 
00262    // -/+inf, which happens to be the Volume default, so we don't have to do 
00263    // anything.
00264    // Otherwise the code to change the bounds would look something like this:
00265 
00266    // first the lower bounds to -inf, upper bounds to inf
00267    volp.dual_lb.allocate(volp.dsize);
00268    volp.dual_lb = -1e31;
00269    volp.dual_ub.allocate(volp.dsize);
00270    volp.dual_ub = 1e31;
00271    // now go through the relaxed constraints and change the lb of the ax >= b 
00272    // constrains to 0, and change the ub of the ax <= b constrains to 0.
00273 /*
00274    for (int i = ufl_data.ncust; i < volp.dsize; ++i) {
00275      volp.dual_ub[i] = 0;  // Assume all (primal) <= ineq.
00276    }
00277 */
00278 
00279 #ifndef WIN32
00280    // start time measurement
00281    double t0;
00282    struct tms timearr; clock_t tres;
00283    tres = times(&timearr); 
00284    t0 = timearr.tms_utime; 
00285 #endif
00286 
00287    // invoke volume algorithm
00288    int vs;
00289    //DEPRECATED
00290    //   if (num_load_balance_masks == 0) {
00291      printf("solving without masks\n");
00292      bool initu = (init_lagrange != 0);
00293      vs = volp.solve(ufl_data, initu);   // JWB ifdual);
00294    /* } else { */
00295 /*      printf("solving with masks\n"); */
00296 /*      vs = volp.solve_with_mask(ufl_data,   */
00297 /*                             load_balance_masks, */
00298 /*                             load_balance_totals, */
00299 /*                             load_balance_parallelize_flags,     */
00300 /*                             num_load_balance_masks, */
00301 /*                             ifdual); */
00302 /*    } */
00303 
00304    if (vs == 0) {
00305       printf("no integer solution during vol run...\n");
00306    } else if (vs == -1) {       // done
00307       printf("solve succeeded within gap...\n");
00308    } else {
00309       // recompute the violation
00310       const int n = ufl_data.nloc;
00311       const int m = ufl_data.ncust;
00312 
00313       VOL_dvector v(volp.dsize);
00314       const VOL_dvector& psol = volp.psol;
00315       v = 1;
00316 
00317 /***  JWB
00318       int i,j,k=n;
00319       for (j = 0; j < n; ++j){
00320         for (i = 0; i < m; ++i) {
00321           v[i] -= psol[k];
00322           ++k;
00323         }
00324       }
00325 ***/
00326 
00327    int k=0;
00328    adjlist& _dist = ufl_data._dist;
00329    //int n_hd = ufl_data.get_n_hd();
00330    //int n_ld = ufl_data.get_n_ld();
00331    //int *high_degree = ufl_data.get_high_degree();
00332    //int *low_degree  = ufl_data.get_low_degree();
00333    int rows = _dist.MatrixRows();
00334    int total = 0;               // debug
00335 
00336    #pragma mta assert parallel
00337    for (int i=0; i < nloc; i++) {
00338      int i_index = _dist.col_index(i);
00339      int next_index = _dist.col_index(i+1);
00340      column_iterator begin_cols = _dist.col_indices_begin(i);
00341      column_iterator end_cols = _dist.col_indices_end(i);
00342      column_iterator ptr2 = begin_cols;
00343      double pcost_incr = 0.0;
00344 
00345         int offs=0;
00346         for (int k=i_index; k<next_index; k++, offs++) {
00347                 column_iterator my_ptr2 = ptr2;
00348                 my_ptr2.set_position(offs);
00349                 int j = *my_ptr2;
00350                 if (ufl_data.ix[n+k] == 1) {
00351                         server[j] = i;
00352                         //printf("ufl finishing: server[%d]: %d (index %d)\n", j, i,n+k);
00353                 }
00354                 v[j] -= psol[n+k];
00355         }
00356    }
00357 
00358 
00359       double vc = 0.0;
00360       for (int i = 0; i < m; ++i) {
00361          double b = ((v[i] < 0) * -v[i]) + ((v[i] >= 0) * v[i]); // MTA
00362          if (b > 0.0)
00363                 vc += b;
00364       }
00365       vc /= m;
00366       printf(" Average violation of final solution: %f\n", vc);
00367 
00368       if (ufl_par.dual_savefile.length() > 0) {
00369         // save dual solution
00370          FILE* file = fopen(ufl_par.dual_savefile.c_str(), "w");
00371          const VOL_dvector& u = volp.dsol;
00372          int n = u.size();
00373          int i;
00374          fclose(file);
00375       }
00376 
00377       // run a couple more heuristics
00378 /*
00379       double heur_val;
00380       for (int i = 0; i < ufl_par.h_iter; ++i) {
00381          heur_val = DBL_MAX;
00382          ufl_data.heuristics(volp, psol, heur_val,0.0);
00383       }
00384 */
00385       // save integer solution
00386       if (ufl_par.int_savefile.length() > 0) {
00387          FILE* file = fopen(ufl_par.int_savefile.c_str(), "w");
00388          const VOL_ivector& x = ufl_data.ix;
00389          const int n = ufl_data.nloc;
00390          const int m = ufl_data.ncust;
00391          int i,j,k=n;
00392          //   SPOT: Need to check both of these, I just guessed:
00393          //         omit the dummy location in this list
00394          //         we need to print out input IDs, not our IDs
00395          fprintf(file, "Best integer solution value: %f\n", 
00396                                 ufl_data.icost);
00397          fprintf(file, "Lower bound: %f\n", volp.value);
00398          fprintf(file, "Sensor locations\n");
00399 
00400          if (origID.size() > 0) {
00401                 for (i = 0,j=0; i < n; ++i) {
00402                         if ( x[i]==1 && (i != newDummyID)){
00403                                 fprintf(file, "%6d ", origID[i]);
00404                                 if (j && (j%10 == 0)) fprintf(file, "\n");
00405                                 j++;
00406                         } 
00407                 }
00408         } else if (open_facilities) {
00409                 #pragma mta assert nodep
00410                 for (i = 0; i < n; ++i) {
00411                         if ( x[i]==1 && (i != newDummyID)){
00412                                 open_facilities[i] = true;
00413                         } else {
00414                                 open_facilities[i] = false;
00415                         }
00416                 }
00417         }
00418          fclose(file);
00419       }
00420    }
00421 
00422    if (frac_solution) {
00423         #pragma mta assert nodep
00424         for (int i=0; i<nloc+nnz; i++) {
00425                 frac_solution[i] = volp.psol[i];
00426         }
00427    }
00428    printf(" ncust: %d\n", ufl_data.ncust);
00429    printf(" Best integer solution value: %f\n", ufl_data.icost);
00430    printf(" Lower bound: %f\n", volp.value);
00431    
00432 #ifndef WIN32
00433    // end time measurement
00434    tres = times(&timearr);
00435    double t = (timearr.tms_utime-t0)/100.;
00436    printf(" Total Time: %f secs\n", t);
00437 #endif
00438 
00439    return 0;
00440 }
00441 
00442 // Read cost objective data or side constraint data.
00443 // Output: Sparse matrix, Matrix.
00444 template <typename vtype, typename adjlist_t>
00445 void UFL_read_data(const char* fname, UFL<vtype,adjlist_t>& data,
00446                    adjlist_t& Matrix) 
00447 {
00448    FILE * file = fopen(fname, "r");
00449    if (!file) {
00450       printf("Failure to open ufl datafile: %s\n ", fname);
00451       abort();
00452    }
00453 
00454 
00455    VOL_dvector& fcost = data.fcost;
00456    //VOL_dvector& dist = data.dist;
00457    //SparseMatrixCSR<double>& _dist = data._dist;
00458 
00459    int& nloc = data.nloc;
00460    int& ncust = data.ncust;
00461    int& nnz = data.nnz;
00462    int& p = data._p;
00463    int len;
00464    int nloc2, ncust2, nnz2;
00465 #if 1
00466    char s[500];
00467    fgets(s, 500, file);
00468    len = strlen(s) - 1;
00469    if (s[len] == '\n')
00470       s[len] = 0;
00471 
00472    // read number of locations, number of customers, and no. of nonzeros
00473    sscanf(s,"%d%d%d",&nloc2,&ncust2,&nnz2);
00474    // check consistency, all constraints must have same sparsity pattern!
00475    // TODO: We only check nnz for now, since nloc may shrink after ID mapping
00476    if (nloc==0)
00477      nloc = nloc2;
00478 #if 0
00479    else if (nloc2 != nloc){
00480      // TODO Standardize error handling 
00481      printf("FATAL: nloc=%i is inconsistent with previous value %i\n",
00482        nloc2, nloc);
00483      exit(-1);
00484    }
00485 #endif
00486    if (ncust==0)
00487      ncust = ncust2;
00488 #if 0
00489    else if (ncust2 != ncust){
00490      // TODO Standardize error handling 
00491      printf("FATAL: ncust=%i is inconsistent with previous value %i\n",
00492        ncust2, ncust);
00493      exit(-1);
00494    }
00495 #endif
00496    if (nnz==0)
00497      nnz = nnz2;
00498    else if (nnz2 != nnz){
00499      // TODO Standardize error handling 
00500      printf("FATAL: nnz=%i is inconsistent with previous value %i\n",
00501        nnz2, nnz);
00502      exit(-1);
00503    }
00504 
00505    fcost.allocate(nloc);
00506    //dist.allocate(nloc*ncust);
00507    int    *m_index  =  (int*)    malloc((nloc+1) * sizeof(int));
00508    double *m_values =  (double*) malloc(nnz    * sizeof(double));
00509    int    *m_columns = (int*)    malloc(nnz    * sizeof(int));
00510    double cost;
00511    int i,j,k;
00512    if (p==0) { // UFL 
00513      // read location costs
00514      for (i = 0; i < nloc; ++i) { 
00515        fgets(s, 500, file);
00516        len = strlen(s) - 1;
00517        if (s[len] == '\n')
00518         s[len] = 0;
00519        sscanf(s,"%lf",&cost);
00520        fcost[i]=cost;
00521      }
00522    }
00523    else { // p-median
00524      // No facility open costs; set them to zero.
00525      for (i = 0; i < nloc; ++i) 
00526        fcost[i]=0.0;
00527    }
00528 
00529    //dist=1.e7;
00530    nnz = 0;
00531    int nextID = 0;
00532    int newi = 0;
00533    char *fgr = fgets(s, 500, file);
00534    k=sscanf(s,"%d%d%lf",&i,&j,&cost);
00535    assert(k==3);
00536    int row = 0;
00537    int cur_loc = i;
00538    int next_loc = i;
00539    int nz_count = 0;
00540    m_index[row] = 0;
00541    while(fgr) {
00542      len = strlen(s) - 1;
00543      if (s[len] == '\n')
00544         s[len] = 0;
00545      if (next_loc != cur_loc) {
00546         row++;
00547         m_index[row] = nz_count;
00548         cur_loc = next_loc;
00549       }
00550      // read cost of serving a customer from a particular location
00551      // i=location, j=customer
00552 
00553      // SPOT:
00554      // Create a set of location IDs that is consecutive beginning at 0.
00555      // Customer (event) IDs are already consecutive beginning at 1.
00556      // Dummy location is "nloc".
00557      // "p" (number of sensors) was incremented by 1 because dummy
00558      // location will be one of the those found by ufl.
00559 
00560      newIDiter = newID.find(i);
00561      if (newIDiter == newID.end()){
00562 
00563        newi = newID[i] = nextID;
00564 
00565        origID.push_back(i);
00566        if (i == nloc){         // sp sends dummy location as nloc
00567          newDummyID = nextID;
00568        }
00569        nextID++;
00570      }
00571      else{
00572        newi = newIDiter->second;
00573      }
00574 
00575      m_values[nz_count] = cost;
00576      m_columns[nz_count] = j-1;
00577 
00578      ++nnz;  
00579      nz_count++;
00580 
00581      fgr = fgets(s, 500, file);
00582      k=sscanf(s,"%d%d%lf",&i,&j,&cost);
00583      if(k!=3) break;
00584      if(i==-1)break;
00585      next_loc = i;
00586    }
00587    row++;
00588    m_index[row] = nz_count;
00589    newID.erase(newID.begin(), newID.end());
00590 
00591    nloc = row;
00592 
00593    Matrix.init(nloc,ncust,nnz,m_index,m_values,m_columns);
00594 
00595 // DEBUG
00596 /*
00597    for (int i=0; i<nloc; i++) {
00598      vtype * begin_vals = Matrix.col_values_begin(i);
00599      int    * begin_cols = Matrix.col_indices_begin(i);
00600      vtype * end_vals = Matrix.col_values_end(i);
00601      int    * end_cols = Matrix.col_indices_end(i);
00602      vtype *ptr1 = begin_vals; 
00603      int    *ptr2 = begin_cols;
00604      for (; ptr1 < end_vals; ptr1++, ptr2++) {
00605        int j = *ptr2;
00606        vtype distij = *ptr1;
00607        printf("(%d:%d) ", j, (int) distij);
00608      }
00609      printf("\n");
00610    }
00611 */
00612 // END DEBUG
00613 
00614 #else // DEAD CODE!
00615    fscanf(file, "%i%i", &ncust, &nloc);
00616    fcost.allocate(nloc);
00617    //dist.allocate(nloc*ncust);
00618    int i,j;
00619    for ( j=0; j<ncust; ++j){
00620      for ( i=0; i<nloc; ++i){
00621        fscanf(file, "%f", &Matrix[i][j]);
00622      }
00623    }
00624    for ( i=0; i<nloc; ++i)
00625      fscanf(file, "%f", &fcost[i]);
00626 #endif
00627    fclose(file);
00628 
00629    data.fix.allocate(nloc);
00630    data.fix = -1;
00631    //if (!verifyMap(dist,ncust,Matrix))  
00632 //      fprintf(stderr, "INCORRECT MAP\n");
00633 }
00634 
00635 
00636 
00637 template <typename vtype, typename adjlist_t>
00638 void UFL<vtype,adjlist_t>::init_u(VOL_dvector& u)
00639 {
00640         printf("UFL::init_u called\n");
00641         if (!init_lagrange) return;
00642         int usize = u.size();
00643         #pragma mta assert nodep
00644         for (int i=0; i<usize; i++) {
00645                 u[i] = init_lagrange[i];
00646         }
00647 }
00648 
00649 template <typename vtype, typename adjlist_t>
00650 double UFL<vtype,adjlist_t>::open_more_servers(VOL_ivector& nsol, 
00651                                 int*unserved, int num_unserved)
00652 {
00653    printf("open_more_servers()\n");
00654    int num_extra=0;
00655    double extra_cost = 0.0;
00656    int *best_server = new int[num_unserved];
00657    memset(best_server, -1, num_unserved*sizeof(int));
00658    #pragma mta assert parallel
00659    for (int i=0; i<nloc; i++) {
00660         int i_index = _dist.col_index(i);
00661         int next_index = _dist.col_index(i+1);
00662         column_iterator begin_cols = _dist.col_indices_begin(i);
00663         column_iterator end_cols = _dist.col_indices_end(i);
00664         column_iterator ptr2 = begin_cols;
00665         double best_cost = fcost[i];
00666         for (int k=i_index; k<next_index; k++) {
00667                 column_iterator my_ptr2 = ptr2;
00668                 my_ptr2.set_position(k-i_index);
00669                 int j = *my_ptr2;
00670                 for (int p=0; p<num_unserved; p++) {
00671 #ifdef DEBUG
00672                         if (j == unserved[p]) {
00673                                 printf("ufl: %d is unserved[%d]\n", 
00674                                                 unserved[p], p);
00675                                 fflush(stdout);
00676                                 printf("ufl: best_server[%d] is %d\n", 
00677                                                 p, best_server[p]);
00678                                 fflush(stdout);
00679                                 if (best_server[p] >= 0) {
00680                                 printf("ufl: fcost[best_server[%d]] is %f\n", 
00681                                                 p, fcost[best_server[p]]);
00682                                 fflush(stdout);
00683                                 printf("ufl: fcost[server %d] is %f\n", 
00684                                                 i, fcost[i]);
00685                                 fflush(stdout);
00686                                 }
00687                         }
00688 #endif
00689                         if ((j == unserved[p]) &&
00690                             ((best_server[p] == -1) ||
00691                              (fcost[i] < fcost[best_server[p]]))) {
00692                                         best_server[p] = i;
00693                         }
00694                 }
00695         }
00696    }
00697    #pragma mta assert parallel
00698    for (int p=0; p<num_unserved; p++) {
00699         int uns = unserved[p];
00700         int serv= best_server[p];
00701         int token = mt_incr(nsol[serv], 1);
00702         if (token == 0) {
00703                 extra_cost += fcost[serv];
00704         } else {
00705                 mt_incr(nsol[serv], -1);
00706         }
00707    }
00708    return extra_cost;
00709 }
00710 
00711 template <typename vtype, typename adjlist_t>
00712 int * UFL<vtype,adjlist_t>::find_unserved(VOL_ivector& nsol, int num_servers,
00713                                           int& num_unserved)
00714 {
00715    int *server = new int[num_servers];
00716    int rows = _dist.MatrixRows();
00717    bool *unserved = new bool[ncust];
00718    int next_server_index = 0;
00719 
00720    #pragma mta assert nodep
00721    for (int i=0; i<ncust; i++) {
00722         unserved[i] = true;
00723    }
00724    #pragma mta assert parallel
00725    for (int i=0; i<rows; i++) {
00726         if (nsol[i] == 1) {
00727                 int next = mt_incr(next_server_index, 1);
00728                 server[next] = i;
00729         }
00730    }
00731    #pragma mta assert parallel
00732    #pragma mta loop future
00733    for (int i=0; i<num_servers; i++) {
00734      int i_index = _dist.col_index(server[i]);
00735      int next_index = _dist.col_index(server[i]+1);
00736      column_iterator begin_cols = _dist.col_indices_begin(server[i]);
00737      column_iterator ptr2 = begin_cols;
00738      //printf("funs: server %d of %d: (%d-%d)\n", i, num_servers,next_index,
00739 //                      i_index);
00740      //fflush(stdout);
00741      if (next_index-i_index > PARCUTOFF) {
00742         #pragma mta assert parallel
00743         for (int k=i_index; k<next_index; k++) {
00744                 int offs = k-i_index;
00745                 ptr2.set_position(offs);
00746                 int j = *ptr2;
00747                 unserved[j] = false;
00748         }
00749      } else {
00750         for (int k=i_index; k<next_index; k++) {
00751                 //printf("funs: server %d of %d: adj %d\n", i, 
00752                 //              num_servers,k);
00753                 fflush(stdout);
00754                 int offs = k-i_index;
00755                 ptr2.set_position(offs);
00756                 int j = *ptr2;
00757                 unserved[j] = false;
00758         }
00759      }
00760    }
00761    int ucount=0;
00762    #pragma mta assert parallel
00763    for (int j=0; j<ncust; j++) {
00764         if (unserved[j]) { 
00765                 ucount += 1;
00766         }
00767    }
00768    num_unserved = ucount;
00769    int *the_unserved = new int[ucount];
00770    ucount = 0;
00771    #pragma mta assert parallel
00772    for (int j=0; j<ncust; j++) {
00773         if (unserved[j]) {   
00774                 int ind = mt_incr(ucount,1);
00775                 the_unserved[ind] = j;
00776                 //printf("unserved: %d\n", j);
00777         }
00778   }
00779   printf("find_unserved: unserved: %d\n", ucount);
00780   delete[] unserved;
00781   delete[] server;
00782   return the_unserved;
00783 }
00784 
00785 template <typename vtype, typename adjlist_t>
00786 double UFL<vtype,adjlist_t>::assign_service(VOL_ivector& nsol, 
00787                                             const VOL_dvector& x)
00788 {
00789    int rows = _dist.MatrixRows();
00790    bool *used = new bool[rows];
00791    double *best_serverness = new double[ncust];
00792    int next_server_index = 0;
00793 
00794    #pragma mta assert nodep
00795    for (int i=0; i<rows; i++) {
00796         used[i] = false;
00797    }
00798    #pragma mta assert nodep
00799    for (int i=0; i<ncust; i++) {
00800         best_serverness[i] = 0.0;
00801    }
00802 
00803    #pragma mta assert parallel
00804    for (int i=0; i<nloc; i++) {
00805      int i_index = _dist.col_index(i);
00806      int next_index = _dist.col_index(i+1);
00807      column_iterator begin_cols = _dist.col_indices_begin(i);
00808      column_iterator end_cols = _dist.col_indices_end(i);
00809      column_iterator ptr2 = begin_cols;
00810      double serverness = nsol[i] * x[i];
00811      //printf("assign_serv: serverness[%d]: %f\n", i, serverness);
00812      //printf("assign_serv: %d:  index[%d,%d]\n", i, i_index, next_index);
00813      for (int k=i_index; k<next_index; k++) {
00814                 double my_serverness = serverness;
00815                 column_iterator my_ptr = ptr2;
00816                 my_ptr.set_position(k-i_index);
00817                 int j = *my_ptr;
00818      //         printf("assign_serv: k: %d, pos: %d, (%d,%d): serverness: %f\n", k, k-i_index, i, j,serverness);
00819                 if (my_serverness > best_serverness[j]) {
00820                         best_serverness[j] = my_serverness;
00821 //                      printf("assign_serv: updating best_serv[%d] to be %d (%f)\n",
00822 //                              j,i,my_serverness);
00823                 }
00824      }
00825      int token = 0;
00826      #pragma mta assert parallel
00827      for (int k=i_index; k<next_index; k++) {
00828                 column_iterator my_ptr = ptr2;
00829                 my_ptr.set_position(k-i_index);
00830                 int j = *my_ptr;
00831                 if (nsol[i] && (x[i] == best_serverness[j])) {
00832                         int got_it = mt_incr(token,1);
00833                         if (got_it==0) {
00834                                 used[i] = true;
00835                         }
00836                 }
00837      }
00838    }
00839    int s_ind = 0;
00840    #pragma mta assert parallel
00841    for (int i=0; i<rows; i++) { 
00842         if ((nsol[i] == 1) && used[i]) {
00843                 int next = mt_incr(s_ind, 1);
00844         } else {
00845                 nsol[i] = 0;
00846         }
00847    }
00848    int num_servers = s_ind;
00849    int *served = new int[ncust];
00850    int *server = new int[num_servers];
00851    s_ind = 0;
00852 
00853    memset(served,0,ncust*sizeof(int));
00854 
00855    #pragma mta assert parallel
00856    for (int i=0; i<rows; i++) { 
00857         if (nsol[i] == 1) {
00858                 int next = mt_incr(next_server_index, 1);
00859                 server[next] = i;
00860         }
00861    }
00862    double cost = 0.0;
00863 
00864    #pragma mta assert parallel
00865    for (int ind=0; ind<next_server_index; ind++) {
00866      int i = server[ind];
00867      int i_index = _dist.col_index(i);
00868      int next_index = _dist.col_index(i+1);
00869      column_iterator begin_cols = _dist.col_indices_begin(i);
00870      column_iterator end_cols = _dist.col_indices_end(i);
00871      value_iterator begin_vals = _dist.col_values_begin(i);
00872      value_iterator end_vals = _dist.col_values_end(i);
00873      column_iterator ptr2 = begin_cols;
00874      value_iterator ptr1 = begin_vals;
00875      double tmp_cost = 0;
00876         for (int k=i_index; k<next_index; k++) {
00877                 column_iterator my_ptr = ptr2;
00878                 my_ptr.set_position(k-i_index);
00879                 int j = *my_ptr;
00880                 int got_it = 0;
00881                 got_it = mt_incr(served[j], 1);  
00882                 //printf("assign_serv: %d trying to serve %d\n", i, j);
00883                 if (got_it == 0) {
00884                 //      printf("assign_serv: %d succeeds in serving %d (index %d)\n",i,j,nloc+k);
00885                         ptr1.set_position(k-i_index);
00886                         double distij = *ptr1;
00887                         tmp_cost += distij;
00888                         nsol[nloc+k] = 1;
00889                 }
00890         }
00891      cost += tmp_cost;
00892    }
00893    delete [] used;
00894    delete [] served;
00895    delete [] server;
00896    delete [] best_serverness;
00897    //delete [] hd;
00898    //delete [] ld;
00899    return cost;
00900 }
00901 
00902 // IN:  fractional primal solution (x),
00903 //      best feasible soln value so far (icost)
00904 // OUT: integral primal solution (ix) if better than best so far
00905 //      and primal value (icost)
00906 // returns -1 if Volume should stop, 0/1 if feasible solution wasn't/was
00907 //         found.
00908 // We use randomized rounding. We look at x[i] as the probability
00909 // of opening facility i.
00910 
00911 template <typename vtype, typename adjlist_t>
00912 int
00913 UFL<vtype,adjlist_t>::heuristics(const VOL_problem& p,
00914                 const VOL_dvector& x, double& new_icost, double lb)
00915 {
00916    VOL_ivector nsol(nloc + nnz);
00917    nsol=0;
00918    int i,j;
00919    double r,value=0;
00920    int state = 12345;
00921    bool int_sol_found=false;
00922    double *xmin = new double[ncust];
00923    int *kmin = new int[ncust];
00924    int fail = 0;
00925    int nopen=0;
00926 
00927   if (fix.size() > 0 ) {
00928         printf("WARNING: heuristics currently ignores fixed locations\n");
00929   }
00930 
00931   nopen=0;
00932   value = 0;
00933 
00934   random_value* rand = (random_value*) malloc(nloc * sizeof(random_value));
00935   rand = rand_fill::generate(nloc, rand);
00936 
00937   if (_p > 0) {  // p_median
00938    for (int i=0; i < nloc; ++i){ // open or close facilities
00939      r=rand[i]/(double)rand_fill::get_max();
00940      // Simple randomized rounding
00941      // TODO: Improve for p-median case (Jon & Cindy)
00942      if (r < x[i]) 
00943          if (nopen < _p){
00944            nsol[i]=1;
00945            ++nopen;
00946          }
00947    }
00948    if (nopen < _p) fail = true;
00949   } else {       // UFL
00950    #pragma mta assert parallel
00951    for (int i=0; i < nloc; ++i){ // open or close facilities
00952      r=rand[i]/(double)rand_fill::get_max();
00953      if (r < x[i]) {
00954            nsol[i]=1;
00955            mt_incr(nopen, 1);
00956            value += fcost[i];
00957      }
00958    }
00959    int num_unserved = 0;
00960    int *the_unserved = find_unserved(nsol, nopen, num_unserved);
00961    value += open_more_servers(nsol, the_unserved, num_unserved);
00962    value += assign_service(nsol, x);
00963   }
00964 
00965    new_icost = value;
00966    if (value < icost) {
00967       icost = value;
00968       ix = nsol;
00969     //printf("heuristics recording a feasible solution:\n");
00970     //for (int i=0; i<nloc; i++) {
00971 //      printf("s[%d] == %d\n", i, ix[i]);
00972     //}
00973    }
00974    /*
00975    printf("ICOST: %d\n", icost);
00976    printf("LB   : %lf\n",lb);
00977    printf("GAP  : %lf\n",(icost-lb)/lb);
00978    */
00979    if (lb > icost) {
00980         printf("lb: %lf exceeds icost: %lf\n", lb, icost);
00981    }
00982    if (lb > 0 && ((icost - lb)/lb < gap) || (lb == 0 && icost < gap)) {
00983         printf("HEURISTICS SAYS TO QUIT! (%lf - %lf)/%lf < %lf\n",
00984                         icost, lb, lb, gap);
00985         return -1;
00986    }
00987    //printf("int sol %f\n", new_icost);
00988    delete [] xmin;
00989    delete [] kmin;
00990 
00991    return 1;
00992 }
00993 
00994 // We assume the side constraint sparsity pattern is the same as 
00995 // the costs _dist, without checking.
00996 // TODO: If not, we need merge the two sparsity patterns.
00997 // This affects the length of the rc vector!
00998 template <typename vtype, typename adjlist_t>
00999 int
01000 UFL<vtype,adjlist_t>::compute_rc(const VOL_dvector& u, VOL_dvector& rc)
01001 {
01002    int i,j,k=0;
01003 
01004    //printf("Debug in comp_rc: Multipliers u (size=%i) = \n", u.size());
01005    //for (i=0; i<u.size(); i++)
01006    //  printf("%lf ", u[i]);
01007    //printf("\n");
01008  
01009 // Standard UFL or p-median objective
01010 
01011    #pragma mta assert parallel
01012    #pragma mta loop future
01013    for ( i=0; i < nloc; i++){
01014      rc[i]= fcost[i];
01015      int i_index = _dist.col_index(i);
01016      int next_index = _dist.col_index(i+1);
01017      value_iterator  begin_vals = _dist.col_values_begin(i);
01018      column_iterator begin_cols = _dist.col_indices_begin(i);
01019      value_iterator  end_vals = _dist.col_values_end(i);
01020      column_iterator end_cols = _dist.col_indices_end(i);
01021      value_iterator  ptr1 = begin_vals; 
01022      column_iterator ptr2 = begin_cols;
01023      if (next_index-i_index > PARCUTOFF) {
01024         int offs = 0;
01025         #pragma mta assert parallel
01026         for (int k=i_index; k<next_index; k++, offs++) {
01027                 value_iterator my_ptr1 = ptr1;
01028                 column_iterator my_ptr2 = ptr2;
01029                 my_ptr1.set_position(offs);
01030                 my_ptr2.set_position(offs);
01031                 int j = *my_ptr2;
01032                 vtype distij = *my_ptr1;
01033                 rc[nloc+k]= distij - u[j];
01034         }
01035      } else {
01036         for (int k=i_index; k<next_index; k++, ptr1++, ptr2++) {
01037                 int j = *ptr2;
01038                 vtype distij = *ptr1;
01039                 rc[nloc+k]= distij - u[j];
01040         }
01041      }
01042    }
01043 
01044    // Additional (linear) side constraints
01045    int c=0;
01046    typename vector<adjlist_t>::iterator sc; 
01047    for (sc = _sidecon.begin(), c=0; sc != _sidecon.end(); ++sc,++c) {
01048      k= 0;
01049      double mult = u[ncust+c]; // Multiplier for side constraint c
01050      //printf("Compute_rc for constraint %d: lambda = %f\n", c, lambda);
01051      for ( i=0; i < nloc; i++){
01052        value_iterator begin_vals = sc->col_values_begin(i);
01053        value_iterator end_vals = sc->col_values_end(i);
01054        value_iterator ptr1 = begin_vals;
01055        for (; ptr1 < end_vals; ptr1++) {
01056          rc[nloc+k] += mult * (*ptr1);
01057          ++k;
01058        }
01059      }
01060    }
01061 
01062    return 0;
01063 }
01064 
01065 /*
01066 template <typename vtype, typename adjlist_t>
01067 void 
01068 UFL<vtype,adjlist_t>::separate_by_degree()
01069 {
01070         int *index = _dist.get_index();
01071         
01072         if (twenty_fifth == 0) {
01073           int *degs = new int[nloc];
01074      
01075           #pragma mta assert nodep
01076           for (int i=0; i<nloc; i++) {
01077             degs[i] = index[i+1] - index[i];
01078           }
01079           int maxval = 0;
01080           #pragma mta assert nodep
01081           for (int i=0; i<nloc; i++) {
01082             if (degs[i] > maxval)
01083               maxval = degs[i];  // compilers removes reduction
01084           }
01085           mtgl::countingSort(degs, nloc, maxval);
01086           twenty_fifth = degs[nloc-25];
01087           //twenty_fifth = 5000;
01088           delete [] degs;
01089           printf("set twenty_fifth to %d vs PARCUTOFF of %d\n", twenty_fifth, PARCUTOFF);
01090         }
01091         
01092 
01093         int next_hd =0;
01094         int next_ld =0;
01095         
01096         #pragma mta assert nodep
01097         for (int i=0; i<nloc; i++) {
01098                 if ((index[i+1] - index[i]) > twenty_fifth) {     
01099                         next_hd++;
01100                 } else {
01101                         next_ld++;
01102                 }
01103         }
01104         n_hd = next_hd;
01105         n_ld = next_ld;
01106         if (n_hd > 0) high_degree = new int[n_hd];
01107         if (n_ld > 0) low_degree  = new int[n_ld];
01108         next_hd = next_ld = 0;
01109         #pragma mta assert parallel
01110         for (int i=0; i<nloc; i++) {
01111                 if ((index[i+1] - index[i]) > twenty_fifth) {     
01112                         int ind = mt_incr(next_hd,1);
01113                         high_degree[ind] = i;
01114                 } else {
01115                         int ind = mt_incr(next_ld,1);
01116                         low_degree[ind] = i;
01117                 }
01118         }
01119         printf("n_hd: %d\n", n_hd);
01120         printf("n_ld: %d\n", n_ld);
01121 }
01122 
01123 template <typename vtype, typename adjlist_t>
01124 void 
01125 UFL<vtype,adjlist_t>::separate_by_degree(int *to_sep, int num_to_sep,
01126                                          int *&hd, int& num_hd,
01127                                          int *&ld, int& num_ld)
01128 {
01129         int next_hd =0;
01130         int next_ld =0;
01131         int *index = _dist.get_index();
01132         #pragma mta assert nodep
01133         for (int ind=0; ind<num_to_sep; ind++) {
01134                 int i = to_sep[ind];
01135                 if ((index[i+1] - index[i]) > PARCUTOFF) {
01136                         next_hd++;
01137                 } else {
01138                         next_ld++;
01139                 }
01140         }
01141         num_hd = next_hd;
01142         num_ld = next_ld;
01143         hd = new int[num_hd];
01144         ld  = new int[num_ld];
01145         next_hd = next_ld = 0;
01146         #pragma mta assert parallel
01147         for (int ind=0; ind<num_to_sep; ind++) {
01148                 int i = to_sep[ind];
01149                 if ((index[i+1] - index[i]) > PARCUTOFF) {
01150                         int ind = mt_incr(next_hd,1);
01151                         hd[ind] = i;
01152                 } else {
01153                         int ind = mt_incr(next_ld,1);
01154                         ld[ind] = i;
01155                 }
01156         }
01157 }
01158 */
01159                 
01160 
01161 // IN: dual vector u
01162 // OUT: primal solution to the Lagrangian subproblem (x)
01163 //      optimal value of Lagrangian subproblem (lcost)
01164 //      v = difference between the rhs and lhs when substituting
01165 //                  x into the relaxed constraints (v)
01166 //      objective value of x substituted into the original problem (pcost)
01167 //      xrc
01168 // return value: -1 (volume should quit) 0 infeasible 1 feasible
01169 
01170 template <typename vtype, typename adjlist_t>
01171 int 
01172 UFL<vtype,adjlist_t>::solve_subproblem(const VOL_dvector& u, 
01173                                        const VOL_dvector& rc,
01174                                        double& lcost, 
01175                                   VOL_dvector& x, VOL_dvector& v, double& pcost)
01176 {
01177    int i,j;
01178 
01179    double lc=0.0;
01180    double *u_v = u.v;
01181    for (i = 0; i < ncust; ++i) {
01182       double incr = u_v[i];
01183       lc += incr;
01184    }
01185    lcost = lc;
01186 
01187    int my_nloc = nloc;
01188 
01189    // VOL_ivector sol(nloc + nloc*ncust);
01190    VOL_ivector sol(x.size());
01191    int *sol_v = sol.v;
01192    double *fcost_v = fcost.v;
01193 
01194    // Produce a primal solution of the relaxed problem
01195    // For p-median, we can only open p facilities
01196  
01197    const double * rdist = rc.v + my_nloc;
01198    int k=0, k1=0;
01199    double value=0.;
01200    int    num_open=0;
01201    int xi;
01202    std::map<int,double>::iterator col_it;
01203    double sum_fcosts = 0.0;
01204 
01205    if (_p == 0) { // UFL
01206      //int first_loop_concur = mta_get_task_counter(RT_CONCURRENCY);
01207      //mt_timer mttimer;
01208      //mttimer.start();
01209      // first treat the high degree nodes in serial
01210      // NOTE: the version that tried load balance masks is in svn 1330
01211 
01212         mt_timer timer;
01213         int issues, memrefs, concur, streams, traps;
01214 
01215      int *fix_v = fix.v;
01216      double *fcost_v = fcost.v;
01217      int fix_size = fix.size();
01218 
01219      #pragma mta assert parallel
01220      for (int i=0; i < my_nloc; ++i ) {
01221        double sum=0.;
01222        double tmp_value=0.;
01223        int xi;
01224        int i_index = _dist.col_index(i);
01225        int next_index = _dist.col_index(i+1);
01226        #pragma mta assert nodep
01227        for (int k=i_index; k<next_index; k++) { 
01228          if ( rdist[k]<0. ) { 
01229            sum += rdist[k];
01230          }
01231          sol_v[my_nloc+k] = 0;
01232        }
01233        if (fix_size > 0) {
01234          if (fix_v[i]==0) xi=0;
01235          else 
01236            if (fix_v[i]==1) xi=1;
01237            else 
01238              if ( fcost_v[i]+sum >= 0. ) xi=0;
01239              else xi=1;
01240        } else {
01241          if ( fcost_v[i]+sum >= 0. ) 
01242            xi=0;
01243          else 
01244            xi=1;
01245        }
01246        sol_v[i]=xi;
01247        double cost_incr = xi * fcost_v[i];
01248        sum_fcosts += cost_incr;
01249        if (xi == 1) {
01250          mt_incr(num_open,1);
01251          
01252          int i_index = _dist.col_index(i);
01253          int next_index = _dist.col_index(i+1);
01254          column_iterator begin_cols = _dist.col_indices_begin(i);
01255          column_iterator end_cols   = _dist.col_indices_end(i);
01256          column_iterator ptr1 = begin_cols;
01257          #pragma mta assert parallel
01258          for (int k=i_index; k<next_index; k++) {
01259            ptr1.set_position(k-i_index);
01260            int j = *ptr1;
01261            if ( rdist[k] < 0. ) {
01262              ptr1.set_position(k-i_index);
01263              int j = *ptr1;
01264              sol_v[my_nloc+k]=(j+1); // remember who's served
01265              tmp_value += rdist[k];
01266              // value += rdist[k] gives internal compiler error with loop future
01267            } 
01268          }
01269        }
01270        value += tmp_value;
01271         //sample_mta_counters(timer, issues, memrefs, concur, streams);
01272         //printf("---------------------------------------------\n");
01273         //printf("solve_sub hd loop: \n");
01274         //printf("---------------------------------------------\n");
01275         //print_mta_counters(timer, nnz, issues,memrefs,concur,streams);
01276      }
01277 
01278 
01279         //printf("---------------------------------------------\n");
01280         //printf("solve_sub ld loop: \n");
01281         //printf("---------------------------------------------\n");
01282         //print_mta_counters(timer, nnz, issues,memrefs,concur,streams);
01283      value += sum_fcosts;
01284    }
01285    else { // _p>0, p-median 
01286      // p lowest rho, as in Avella et al.
01287      // (index, value) pairs
01288      std::vector<std::pair<int,double> > rho(_p+1); 
01289      std::pair<int,double> rho_temp; 
01290      int q;
01291 
01292      for (q=0; q<_p; q++) {
01293        rho[q].first = -1;
01294        rho[q].second = 0;
01295      }
01296 
01297      for ( i=0; i < my_nloc; ++i ) {
01298        double sum=0.;
01299        column_iterator begin_cols = _dist.col_indices_begin(i);
01300        column_iterator end_cols   = _dist.col_indices_end(i);
01301        for (column_iterator ptr1 = begin_cols; ptr1 < end_cols; ptr1++) {
01302          if ( rdist[k]<0. ) sum+=rdist[k];
01303          ++k;
01304        }
01305 
01306        // Save p lowest rho values
01307        rho_temp.first = i;
01308        // rho_temp.second = sum - u[i]; // Avelli et al. use y[j] for x[j,j]
01309        rho_temp.second = sum;
01310        rho[_p] = rho_temp;
01311        for (q=_p; q>0; q--){
01312          if (rho[q-1].first == -1 || (rho[q].second < rho[q-1].second)){
01313            rho_temp = rho[q-1];
01314            rho[q-1] = rho[q];
01315            rho[q] = rho_temp;
01316          }
01317          else
01318            break;
01319        }
01320      }
01321 
01322      // Compute x and function value based on smallest rho values
01323      value = 0;
01324      for ( i=0; i < my_nloc+nnz; ++i ) {
01325        sol[i] = 0;
01326      }
01327      // Loop over open locations
01328      int ii;
01329      for ( ii=0; ii < _p; ++ii ) {
01330        i = rho[ii].first;
01331        // Open a facility at location i
01332        sol[i]=1;
01333        //printf("solvesp: Debug: Open facility at location %d, rho= %f\n",
01334        //   i, rho[ii].second);
01335        value+= rho[ii].second;
01336        int i_index = _dist.col_index(i);
01337        int next_index = _dist.col_index(i+1);
01338        column_iterator begin_cols = _dist.col_indices_begin(i);
01339        column_iterator end_cols   = _dist.col_indices_end(i);
01340        column_iterator ptr1 = begin_cols;
01341        for (int k=i_index; k<next_index; k++, ptr1++) {
01342          // Set customer-facility variables
01343          if ( rdist[k-my_nloc] < 0. ) {
01344                 int j = *ptr1;
01345                 //printf("assigning customer %d to %d\n", j, i);
01346                 sol[k] = 1;
01347          }
01348        }
01349      }
01350    }
01351 
01352    lcost += value; 
01353 
01354    double pcost_tmp = 0.0;
01355    x = 0.0;
01356    double *x_v = x.v;
01357    #pragma mta assert parallel
01358    for (i = 0; i < my_nloc; ++i) {
01359      if (_p == 0) { // Include facility open costs for UFL, not p-median
01360        double incr = fcost_v[i] * sol_v[i];
01361        pcost_tmp += incr;
01362      }
01363      x_v[i] = sol_v[i];
01364    }
01365 
01366    // Compute x, v, pcost from sol
01367    k = 0;
01368    #pragma mta assert parallel
01369    for (int i=0; i < my_nloc; i++) {
01370      int i_index = _dist.col_index(i);
01371      int next_index = _dist.col_index(i+1);
01372      value_iterator  begin_vals = _dist.col_values_begin(i);
01373      column_iterator begin_cols = _dist.col_indices_begin(i);
01374      value_iterator  end_vals = _dist.col_values_end(i);
01375      column_iterator end_cols = _dist.col_indices_end(i);
01376      value_iterator  ptr1 = begin_vals; 
01377      column_iterator ptr2 = begin_cols;
01378      double pcost_incr = 0.0;
01379 
01380         int offs=0;
01381         for (int k=i_index; k<next_index; k++, offs++) {
01382                 value_iterator  my_ptr1 = ptr1;
01383                 column_iterator my_ptr2 = ptr2;
01384                 my_ptr1.set_position(offs);
01385                 my_ptr2.set_position(offs);
01386                 int j = *my_ptr2;
01387                 vtype distij = *my_ptr1;
01388                 double incr = distij*sol_v[my_nloc+k];
01389                 double solnk = (sol_v[my_nloc+k] != 0);
01390                 pcost_incr += incr;
01391                 x_v[my_nloc+k]=solnk;
01392         }
01393     pcost_tmp += pcost_incr;
01394    }
01395    pcost = pcost_tmp;
01396 
01397 
01398    int *tmp_v = new int[ncust];
01399    #pragma mta assert nodep
01400    for (int j=0; j<ncust; j++) {
01401         tmp_v[j] = 1;
01402    }
01403 
01404    int upper = my_nloc+nnz;
01405    #pragma mta assert nodep
01406    for (int k=my_nloc; k<upper; k++) {
01407         if (sol_v[k]) { 
01408                 int& v = tmp_v[sol_v[k]-1];
01409                 v -= 1;
01410         }
01411    }
01412 
01413    #pragma mta assert parallel
01414    for (int j=0; j<ncust; j++) {
01415         v[j] = tmp_v[j];
01416    }
01417    delete [] tmp_v;
01418 
01419    // Compute slack/violations of side constraints
01420    #pragma mta assert nodep
01421    for (i=0; i<_ub.size(); i++)
01422      v[ncust+i] = _ub[i]; 
01423 
01424    int c;
01425    typename vector<adjlist_t>::iterator sc;
01426    for (sc= _sidecon.begin(), c=0; sc != _sidecon.end(); ++sc,++c) {
01427      k = 0;
01428      for ( i=0; i < my_nloc; i++){
01429        int i_index = _dist.col_index(i);
01430        int next_index = _dist.col_index(i+1);
01431        value_iterator  begin_vals = sc->col_values_begin(i);
01432        column_iterator begin_cols = sc->col_indices_begin(i);
01433        value_iterator  end_vals = sc->col_values_end(i);
01434        column_iterator end_cols = sc->col_indices_end(i);
01435        value_iterator  ptr1 = begin_vals; 
01436        column_iterator ptr2 = begin_cols;
01437        for (int k=i_index; k<next_index; k++, ptr1++, ptr2++) {
01438          ptr1.set_position(k);
01439          ptr2.set_position(k);
01440          int j = *ptr2;
01441          double distij = *ptr1;
01442          v[ncust+c] -= distij * x[my_nloc+k]; // Update side constraint c
01443        }
01444      }
01445      // Update lcost (Lagrangian value)
01446      lcost -= u[ncust+c] * _ub[c]; // TODO Check + or - ?
01447    }
01448 
01449    return 1;
01450 }
01451 
01452 
01453 
01454 
01455 
01456 // DEPRECATED
01457 // THIS USES THE MASKS
01458 // IN: dual vector u
01459 // OUT: primal solution to the Lagrangian subproblem (x)
01460 //      optimal value of Lagrangian subproblem (lcost)
01461 //      v = difference between the rhs and lhs when substituting
01462 //                  x into the relaxed constraints (v)
01463 //      objective value of x substituted into the original problem (pcost)
01464 //      xrc
01465 // return value: -1 (volume should quit) 0 infeasible 1 feasible
01466 
01467 /* template <typename vtype, typename adjlist_t> */
01468 /* int  */
01469 /* UFL<vtype,adjlist_t>::solve_subproblem_with_mask(const VOL_dvector& u,  */
01470 /*                                     const VOL_dvector& rc, */
01471 /*                                     double& lcost,  */
01472 /*                                     VOL_dvector& x, VOL_dvector& v, double& pcost,       */
01473 /*                                     int **load_balance_masks, */
01474 /*                                     int *load_balance_totals, */
01475 /*                                             int **load_balance_parallelize_flags,     */
01476 /*                                     int num_load_balance_masks */
01477 /*                                               ) */
01478 /* { */
01479 /*    int i,j; */
01480 
01481 /*    double lc=0.0; */
01482 /*    for (i = 0; i < ncust; ++i) { */
01483 /*       double incr = u[i]; */
01484 /*       lc += incr; */
01485 /*    } */
01486 /*    lcost = lc; */
01487 
01488 /*    // VOL_ivector sol(nloc + nloc*ncust); */
01489 /*    VOL_ivector sol(x.size()); */
01490 
01491 /*    // Produce a primal solution of the relaxed problem */
01492 /*    // For p-median, we can only open p facilities */
01493  
01494 /*    const double * rdist = rc.v + nloc; */
01495 /*    int k=0, k1=0; */
01496 /*    double value=0.; */
01497 /*    int    num_open=0; */
01498 /*    //   int xi; */
01499 /*    std::map<int,double>::iterator col_it; */
01500 /*    double sum_fcosts = 0.0; */
01501 
01502 /*    if (_p == 0) { // UFL */
01503 /*      //int first_loop_concur = mta_get_task_counter(RT_CONCURRENCY); */
01504 /*      //nt_timer mttimer; */
01505 /*      //mttimer.start(); */
01506 /*      // first treat the high degree nodes in serial */
01507 /*      // NOTE: the version that tried load balance masks is in svn 1330 */
01508 
01509 /*      //start new code */
01510 
01511 /*      //could there be time savings if we index into the part of sol that we need for each vertex? */
01512 /*      //eg. ptr = sol[offset+k]? */
01513 /*      //also check for hotspots!!! */
01514 /*      //we should also have a better way to pass in the mask info (struct maybe?) */
01515 /*      for (int i = 0; i < num_load_balance_masks; i++) { */
01516 /*        if (load_balance_parallelize_flags[i][0] == SERIAL) {     */
01517 /*       for (int j = 0; j < load_balance_totals[i]; j++){ */
01518 /*         //using tmp variables to see if that would cut down on hotspots */
01519 /*         double tmp_sum_fcosts = 0, tmp_value = 0; */
01520 /*         process_opening_serving_of_nodes(i, j, num_open, tmp_sum_fcosts, tmp_value, sol, rdist,  */
01521 /*                                          load_balance_masks, load_balance_totals,  */
01522 /*                                          load_balance_parallelize_flags, num_load_balance_masks); */
01523 /*         sum_fcosts += tmp_sum_fcosts; */
01524 /*         value      += tmp_value; */
01525 /*       } */
01526 /*        } else if (load_balance_parallelize_flags[i][0] == PARALLEL) {     */
01527 /*        #pragma mta assert parallel */
01528 /*       for (int j = 0; j < load_balance_totals[i]; j++){ */
01529 /*         //using tmp variables to see if that would cut down on hotspots */
01530 /*         double tmp_sum_fcosts = 0, tmp_value = 0; */
01531 /*         process_opening_serving_of_nodes(i, j, num_open, tmp_sum_fcosts, tmp_value, sol, rdist, */
01532 /*                                          load_balance_masks, load_balance_totals,  */
01533 /*                                          load_balance_parallelize_flags, num_load_balance_masks); */
01534 /*         sum_fcosts += tmp_sum_fcosts; */
01535 /*         value      += tmp_value; */
01536 /*       } */
01537 /*        } else if (load_balance_parallelize_flags[i][0] == PARALLEL_FUTURE) {     */
01538 /*          #pragma mta loop future */
01539 /*       for (int j = 0; j < load_balance_totals[i]; j++){ */
01540 /*         //using tmp variables to see if that would cut down on hotspots */
01541 /*         double tmp_sum_fcosts = 0, tmp_value = 0; */
01542 /*         process_opening_serving_of_nodes(i, j, num_open, tmp_sum_fcosts, tmp_value, sol, rdist, */
01543 /*                                          load_balance_masks, load_balance_totals,  */
01544 /*                                          load_balance_parallelize_flags, num_load_balance_masks); */
01545 /*         sum_fcosts += tmp_sum_fcosts; */
01546 /*         value      += tmp_value; */
01547 /*       } */
01548 /*        } */
01549 /*      } */
01550      
01551 /*      value += sum_fcosts; */
01552      
01553 /*      //end new code */
01554 
01555 /*      /\*     for (int ind=0; ind < n_hd; ++ind ) { */
01556 /*        int i = high_degree[ind]; */
01557 /*        k1 = k; */
01558 /*        double sum=0.; */
01559 /*        double tmp_value=0.; */
01560 /*        int i_index = _dist.col_index(i); */
01561 /*        int next_index = _dist.col_index(i+1); */
01562 /*        #pragma mta assert nodep */
01563 /*        for (int k=i_index; k<next_index; k++) {       */
01564 /*          if ( rdist[k]<0. ) {  */
01565 /*         sum += rdist[k]; */
01566 /*       } */
01567 /*       sol[nloc+k] = 0; */
01568 /*        } */
01569 /*        if (fix.size() > 0) { */
01570 /*       if (fix[i]==0) xi=0; */
01571 /*       else  */
01572 /*         if (fix[i]==1) xi=1; */
01573 /*         else  */
01574 /*           if ( fcost[i]+sum >= 0. ) xi=0; */
01575 /*           else xi=1; */
01576 /*        } else { */
01577 /*       if ( fcost[i]+sum >= 0. )  */
01578 /*         xi=0; */
01579 /*       else  */
01580 /*         xi=1; */
01581 /*        } */
01582 /*        sol[i]=xi; */
01583 /*        double cost_incr = xi * fcost[i]; */
01584 /*        sum_fcosts += cost_incr; */
01585 /*        if (xi == 1) { */
01586 /*       mt_incr(num_open,1); */
01587          
01588 /*       int i_index = _dist.col_index(i); */
01589 /*       int next_index = _dist.col_index(i+1); */
01590 /*       column_iterator begin_cols = _dist.col_indices_begin(i); */
01591 /*       column_iterator end_cols   = _dist.col_indices_end(i); */
01592 /*       column_iterator ptr1 = begin_cols; */
01593 /*          #pragma mta assert parallel */
01594 /*       for (int k=i_index; k<next_index; k++) { */
01595 /*         if ( rdist[k] < 0. ) { */
01596 /*           ptr1.set_position(k-i_index); */
01597 /*           int j = *ptr1; */
01598 /*           sol[nloc+k]=(j+1); // remember who's served */
01599 /*           tmp_value += rdist[k]; */
01600 /*           // value += rdist[k] gives internal compiler error with loop future */
01601 /*         }  */
01602 /*       } */
01603 /*        } */
01604 /*        value += tmp_value; */
01605 /*      } */
01606 /*      //int total_first_loop_concur = mta_get_task_counter(RT_CONCURRENCY) - */
01607 /*      //                                      first_loop_concur; */
01608 /*      //mttimer.stop(); */
01609 /*      //int ticks = mttimer.getElapsedTicks(); */
01610 /*      //printf("solve_sub: first loop concur: %lf\n",  */
01611 /*      //                      total_first_loop_concur / (double) ticks); */
01612 /*      #pragma mta assert parallel */
01613 /*      for (int ind=0; ind < n_ld; ++ind ) { */
01614 /*        int i = low_degree[ind]; */
01615 /*        k1 = k; */
01616 /*        double sum=0.; */
01617 /*        double tmp_value=0.; */
01618 /*        int i_index = _dist.col_index(i); */
01619 /*        int next_index = _dist.col_index(i+1); */
01620 /*        #pragma mta assert nodep */
01621 /*        for (int k=i_index; k<next_index; k++) {       */
01622 /*          if ( rdist[k]<0. ) {  */
01623 /*         sum += rdist[k]; */
01624 /*       } */
01625 /*       sol[nloc+k] = 0; */
01626 /*        } */
01627 /*        if (fix.size() > 0) { */
01628 /*       if (fix[i]==0) xi=0; */
01629 /*       else  */
01630 /*         if (fix[i]==1) xi=1; */
01631 /*         else  */
01632 /*           if ( fcost[i]+sum >= 0. ) xi=0; */
01633 /*           else xi=1; */
01634 /*        } else { */
01635 /*       if ( fcost[i]+sum >= 0. )  */
01636 /*         xi=0; */
01637 /*       else  */
01638 /*         xi=1; */
01639 /*        } */
01640 /*        sol[i]=xi; */
01641 /*        double cost_incr = xi * fcost[i]; */
01642 /*        sum_fcosts += cost_incr; */
01643 /*        if (xi == 1) { */
01644 /*       mt_incr(num_open,1); */
01645          
01646 /*       int i_index = _dist.col_index(i); */
01647 /*       int next_index = _dist.col_index(i+1); */
01648 /*       column_iterator begin_cols = _dist.col_indices_begin(i); */
01649 /*       column_iterator end_cols   = _dist.col_indices_end(i); */
01650 /*       column_iterator ptr1 = begin_cols; */
01651 /*       for (int k=i_index; k<next_index; k++) { */
01652 /*         if ( rdist[k] < 0. ) { */
01653 /*           ptr1.set_position(k-i_index); */
01654 /*           int j = *ptr1; */
01655 /*           sol[nloc+k]=(j+1); // remember who's served */
01656 /*           tmp_value += rdist[k]; */
01657 /*           // value += rdist[k] gives internal compiler error with loop future */
01658 /*         }  */
01659 /*       } */
01660 /*        } */
01661 /*        value += tmp_value; */
01662 /*      } */
01663 /*      value += sum_fcosts; */
01664 
01665 /* *\/ */
01666 /*    } */
01667 /*    else { // _p>0, p-median  */
01668 /*      // p lowest rho, as in Avella et al. */
01669 /*      // (index, value) pairs */
01670 /*      std::vector<std::pair<int,double> > rho(_p+1);  */
01671 /*      std::pair<int,double> rho_temp;  */
01672 /*      int q; */
01673 
01674 /*      for (q=0; q<_p; q++) { */
01675 /*        rho[q].first = -1; */
01676 /*        rho[q].second = 0; */
01677 /*      } */
01678 
01679 /*      for ( i=0; i < nloc; ++i ) { */
01680 /*        double sum=0.; */
01681 /*        column_iterator begin_cols = _dist.col_indices_begin(i); */
01682 /*        column_iterator end_cols   = _dist.col_indices_end(i); */
01683 /*        for (column_iterator ptr1 = begin_cols; ptr1 < end_cols; ptr1++) { */
01684 /*          if ( rdist[k]<0. ) sum+=rdist[k]; */
01685 /*          ++k; */
01686 /*        } */
01687 
01688 /*        // Save p lowest rho values */
01689 /*        rho_temp.first = i; */
01690 /*        // rho_temp.second = sum - u[i]; // Avelli et al. use y[j] for x[j,j] */
01691 /*        rho_temp.second = sum; */
01692 /*        rho[_p] = rho_temp; */
01693 /*        for (q=_p; q>0; q--){ */
01694 /*          if (rho[q-1].first == -1 || (rho[q].second < rho[q-1].second)){ */
01695 /*            rho_temp = rho[q-1]; */
01696 /*            rho[q-1] = rho[q]; */
01697 /*            rho[q] = rho_temp; */
01698 /*          } */
01699 /*          else */
01700 /*            break; */
01701 /*        } */
01702 /*      } */
01703 
01704 /*      // Compute x and function value based on smallest rho values */
01705 /*      value = 0; */
01706 /*      for ( i=0; i < nloc+nnz; ++i ) { */
01707 /*        sol[i] = 0; */
01708 /*      } */
01709 /*      // Loop over open locations */
01710 /*      int ii; */
01711 /*      for ( ii=0; ii < _p; ++ii ) { */
01712 /*        i = rho[ii].first; */
01713 /*        // Open a facility at location i */
01714 /*        sol[i]=1; */
01715 /*        //printf("solvesp: Debug: Open facility at location %d, rho= %f\n", */
01716 /*        //   i, rho[ii].second); */
01717 /*        value+= rho[ii].second; */
01718 /*        int i_index = _dist.col_index(i); */
01719 /*        int next_index = _dist.col_index(i+1); */
01720 /*        column_iterator begin_cols = _dist.col_indices_begin(i); */
01721 /*        column_iterator end_cols   = _dist.col_indices_end(i); */
01722 /*        column_iterator ptr1 = begin_cols; */
01723 /*        for (int k=i_index; k<next_index; k++, ptr1++) { */
01724 /*          // Set customer-facility variables */
01725 /*          if ( rdist[k-nloc] < 0. ) { */
01726 /*              int j = *ptr1; */
01727 /*              //printf("assigning customer %d to %d\n", j, i); */
01728 /*              sol[k] = 1; */
01729 /*       } */
01730 /*        } */
01731 /*      } */
01732 /*    } */
01733 
01734 /*    lcost += value;  */
01735 
01736 /*    double pcost_tmp = 0.0; */
01737 /*    x = 0.0; */
01738 /*    #pragma mta assert parallel */
01739 /*    for (i = 0; i < nloc; ++i) { */
01740 /*      if (_p == 0) { // Include facility open costs for UFL, not p-median */
01741 /*        double incr = fcost[i] * sol[i]; */
01742 /*        pcost_tmp += incr; */
01743 /*      } */
01744 /*      x[i] = sol[i]; */
01745 /*    } */
01746 
01747 /*    // Compute x, v, pcost from sol */
01748 /*    k = 0; */
01749 /*    #pragma mta assert parallel */
01750 /*    for (int i=0; i < nloc; i++) { */
01751 /*      int i_index = _dist.col_index(i); */
01752 /*      int next_index = _dist.col_index(i+1); */
01753 /*      value_iterator  begin_vals = _dist.col_values_begin(i); */
01754 /*      column_iterator begin_cols = _dist.col_indices_begin(i); */
01755 /*      value_iterator  end_vals = _dist.col_values_end(i); */
01756 /*      column_iterator end_cols = _dist.col_indices_end(i); */
01757 /*      value_iterator  ptr1 = begin_vals;  */
01758 /*      column_iterator ptr2 = begin_cols; */
01759 /*      double pcost_incr = 0.0; */
01760 
01761 /*      int offs=0; */
01762 /*              for (int k=i_index; k<next_index; k++, offs++) { */
01763 /*              ptr1.set_position(offs); */
01764 /*              ptr2.set_position(offs); */
01765 /*                      int j = *ptr2; */
01766 /*                      vtype distij = *ptr1; */
01767 /*              double incr = distij*sol[nloc+k]; */
01768 /*              double solnk = (sol[nloc+k] != 0); */
01769 /*                      pcost_incr += incr; */
01770 /*                      x[nloc+k]=solnk; */
01771 /*              } */
01772 /*     pcost_tmp += pcost_incr; */
01773 /*    } */
01774 /*    pcost = pcost_tmp; */
01775 
01776 
01777 /*    int *tmp_v = new int[ncust]; */
01778 /*    #pragma mta assert nodep */
01779 /*    for (int j=0; j<ncust; j++) { */
01780 /*      tmp_v[j] = 1; */
01781 /*    } */
01782 
01783 /*    int upper = nloc+nnz; */
01784 /*    #pragma mta assert nodep */
01785 /*    for (int k=nloc; k<upper; k++) { */
01786 /*      if (sol[k]) {  */
01787 /*              int& v = tmp_v[sol[k]-1]; */
01788 /*              v -= 1; */
01789 /*      } */
01790 /*    } */
01791 
01792 /*    #pragma mta assert parallel */
01793 /*    for (int j=0; j<ncust; j++) { */
01794 /*      v[j] = tmp_v[j]; */
01795 /*    } */
01796 /*    delete [] tmp_v; */
01797 
01798 /*    // Compute slack/violations of side constraints */
01799 /*    #pragma mta assert nodep */
01800 /*    for (i=0; i<_ub.size(); i++) */
01801 /*      v[ncust+i] = _ub[i];  */
01802 
01803 /*    int c; */
01804 /*    typename vector<adjlist_t>::iterator sc; */
01805 /*    for (sc= _sidecon.begin(), c=0; sc != _sidecon.end(); ++sc,++c) { */
01806 /*      k = 0; */
01807 /*      for ( i=0; i < nloc; i++){ */
01808 /*        int i_index = _dist.col_index(i); */
01809 /*        int next_index = _dist.col_index(i+1); */
01810 /*        value_iterator  begin_vals = sc->col_values_begin(i); */
01811 /*        column_iterator begin_cols = sc->col_indices_begin(i); */
01812 /*        value_iterator  end_vals = sc->col_values_end(i); */
01813 /*        column_iterator end_cols = sc->col_indices_end(i); */
01814 /*        value_iterator  ptr1 = begin_vals;  */
01815 /*        column_iterator ptr2 = begin_cols; */
01816 /*        for (int k=i_index; k<next_index; k++, ptr1++, ptr2++) { */
01817 /*       ptr1.set_position(k); */
01818 /*       ptr2.set_position(k); */
01819 /*          int j = *ptr2; */
01820 /*          double distij = *ptr1; */
01821 /*          v[ncust+c] -= distij * x[nloc+k]; // Update side constraint c */
01822 /*        } */
01823 /*      } */
01824 /*      // Update lcost (Lagrangian value) */
01825 /*      lcost -= u[ncust+c] * _ub[c]; // TODO Check + or - ? */
01826 /*    } */
01827 
01828 /*    return 1; */
01829 /* } */
01830 
01831 template <typename vtype, typename adjlist_t>
01832   void UFL<vtype,adjlist_t>::process_opening_serving_of_nodes(int mask_iteration, int mask_vertex_index, int &num_open, double &sum_fcosts, double &value, VOL_ivector &sol, const double *rdist, int **load_balance_masks, int *load_balance_totals, int **load_balance_parallelize_flags, int num_load_balance_masks) {
01833   
01834   int vertex_id = load_balance_masks[mask_iteration][mask_vertex_index];
01835   int i_index = _dist.col_index(vertex_id);
01836   int next_index = _dist.col_index(vertex_id + 1);
01837   int xi;
01838   double sum = 0.0;
01839 
01840   if (load_balance_parallelize_flags[mask_iteration][1] == SERIAL) {    
01841     for (int k = i_index; k < next_index; k++){
01842       if (rdist[k] < 0.0) {
01843         sum += rdist[k];
01844       }
01845       sol[nloc+k] = 0;
01846     }
01847   } else if (load_balance_parallelize_flags[mask_iteration][1] == PARALLEL) {    
01848   #pragma mta assert parallel
01849     for (int k = i_index; k < next_index; k++){   
01850       if (rdist[k] < 0.0) {
01851         sum += rdist[k];
01852       }
01853       sol[nloc+k] = 0;
01854     }
01855   } else if (load_balance_parallelize_flags[mask_iteration][1] == PARALLEL_FUTURE) {    
01856   #pragma mta loop future
01857     for (int k = i_index; k < next_index; k++){
01858       if (rdist[k] < 0.0) {
01859         sum += rdist[k];
01860       }
01861       sol[nloc+k] = 0;  
01862     }
01863   }
01864 
01865   if (fix.size() > 0) {
01866     if (fix[vertex_id] == 0)
01867       xi = 0;
01868     else if (fix[vertex_id] == 1)
01869       xi = 1;
01870     else if (fcost[vertex_id] + sum >= 0.0)
01871       xi = 0;
01872     else
01873       xi = 1;
01874   } else {
01875     if (fcost[vertex_id] + sum >= 0.0)
01876       xi = 0;
01877     else
01878       xi = 1;
01879   }
01880 
01881   sol[vertex_id] = xi;
01882   double cost_incr = xi * fcost[vertex_id];
01883 
01884   //TODO: make sure this is a safe addition
01885   sum_fcosts += cost_incr;
01886 
01887   //now update neighbors
01888   double tmp_value = 0.0;
01889   if (xi == 1) {
01890     mt_incr(num_open, 1);
01891     int i_index = _dist.col_index(vertex_id);
01892     int next_index = _dist.col_index(vertex_id+1);
01893     column_iterator ptr1 = _dist.col_indices_begin(vertex_id);
01894 
01895 
01896     // TODO: make sure that these are multithreaded safe additions
01897     if (load_balance_parallelize_flags[mask_iteration][1] == SERIAL) {    
01898       for (int k = i_index; k < next_index; k++){
01899         if (rdist[k] < 0.0) {
01900           ptr1.set_position(k - i_index);
01901           int j = *ptr1;
01902           sol[nloc+k] = (j+1);
01903           tmp_value += rdist[k];
01904         }
01905       }
01906     } else if (load_balance_parallelize_flags[mask_iteration][1] == PARALLEL) {    
01907       #pragma mta assert parallel
01908       for (int k = i_index; k < next_index; k++){   
01909         if (rdist[k] < 0.0) {
01910           ptr1.set_position(k - i_index);
01911           int j = *ptr1;
01912           sol[nloc+k] = (j+1);
01913           tmp_value += rdist[k];
01914         }
01915       }
01916     } else if (load_balance_parallelize_flags[mask_iteration][1] == PARALLEL_FUTURE) {    
01917       #pragma mta loop future
01918       for (int k = i_index; k < next_index; k++){
01919         if (rdist[k] < 0.0) {
01920           ptr1.set_position(k - i_index);
01921           int j = *ptr1;
01922           sol[nloc+k] = (j+1);
01923           tmp_value += rdist[k];
01924         }
01925       }
01926     }
01927     
01928     value += tmp_value;
01929 
01930   }
01931   
01932 }
01933 
01934 #endif

Generated on Fri Oct 22 2010 11:02:23 for SST by  doxygen 1.7.1