Ejemplo n.º 1
0
void Transformer::assignHandler_(CGraphPtr cg, ConstraintPtr c)
{
  switch (cg->getOut()->getOp()) {
  case OpMult:
  case OpSqr:
    qHandler_->addConstraint(c);
    break;
  default:
    {
    VariablePtr iv;
    VariablePtr ov = VariablePtr();
    LinearFunctionPtr lf = c->getFunction()->getLinearFunction();
    if (lf) {
      assert(lf->getNumTerms()==1);
      ov = lf->termsBegin()->first;
    }
    iv = *(c->getFunction()->getNonlinearFunction()->varsBegin());
    uHandler_->addConstraint(c, iv, ov, 'E');
    }
  }
}
Ejemplo n.º 2
0
void ParQGHandler::addCut_(const double *nlpx, const double *lpx,
                        ConstraintPtr con, CutManager *cutman,
                        SeparationStatus *status)
{
  int error=0;
  ConstraintPtr newcon;
  std::stringstream sstm;
  double c, lpvio, act, cUb;
  FunctionPtr f = con->getFunction();
  LinearFunctionPtr lf = LinearFunctionPtr();

  act = con->getActivity(nlpx, &error);
  if (error == 0) {
    linearAt_(f, act, nlpx, &c, &lf, &error);
    if (error==0) {
      cUb = con->getUb();
      lpvio = std::max(lf->eval(lpx)-cUb+c, 0.0);
      if ((lpvio>solAbsTol_) &&
          ((cUb-c)==0 || (lpvio>fabs(cUb-c)*solRelTol_))) {
#if SPEW
        logger_->msgStream(LogDebug) << me_ << " linearization of constraint "
          << con->getName() << " violated at LP solution with violation = " <<
          lpvio << ". OA cut added." << std::endl;
#endif
        ++(stats_->cuts);
        sstm << "_OAcut_";
        sstm << stats_->cuts;
        *status = SepaResolve;
        f = (FunctionPtr) new Function(lf);
        newcon = rel_->newConstraint(f, -INFINITY, cUb-c, sstm.str());
        CutPtr cut = (CutPtr) new Cut(minlp_->getNumVars(),f, -INFINITY,
                                      cUb-c, false,false);
        cut->setCons(newcon);
        cutman->addCutToPool(cut);
        return;
      }
    }
  }	else {
    logger_->msgStream(LogError) << me_ << "Constraint not defined at"
      << " this point. "<<  std::endl;
#if SPEW
          logger_->msgStream(LogDebug) << me_ << "constraint " <<
            con->getName() << " not defined at this point." << std::endl;
#endif
  }
  return;
}
Ejemplo n.º 3
0
void
MultilinearFormulation::makeConvexHull()
{
  
  // First we do some processing of the instance to determine how many 
  //  multilinear or quadratic constraints, and we store the indicies
  // in the instance for later    
  vector<int> lcid;
  vector<int> mlcid;
  for(ConstConstraintIterator it = originalInstance_->consBegin(); 
      it != originalInstance_->consEnd(); ++it) {
    FunctionType ft = (*it)->getFunctionType();
    if (ft == Multilinear) {
      mlcid.push_back((*it)->getId());
    }    
    else if (ft == Bilinear) {
      mlcid.push_back((*it)->getId());
    }
    else if (ft == Linear) {
      lcid.push_back((*it)->getId());
    }
  }
    

  // add x variables
  vector<double> lb;
  vector<double> ub;

  vector<VariablePtr> xvars;
  int nv = 0;
  for (ConstVariableIterator it = originalInstance_->varsBegin(); 
       it != originalInstance_->varsEnd(); ++it) {
    VariablePtr v = *it;
    VariablePtr vnew = VariablePtr(new Variable(nv, v->getLb(), v->getUb(), v->getType()));
    lb.push_back(v->getLb());
    ub.push_back(v->getUb());

    variableMapping_.insert(make_pair(vnew,v));    
    variables_.push_back(vnew);
    xvars.push_back(vnew);
    nv++;
  }

  // Add z vars
  vector<VariablePtr> zvars;
  for(UInt i = 0; i < mlcid.size(); i++) {
    const ConstraintPtr mlc = originalInstance_->getConstraint(mlcid[i]);

    VariablePtr vnew = VariablePtr(new Variable(nv,mlc->getLb(),mlc->getUb(),Continuous));
    variables_.push_back(vnew);
    zvars.push_back(vnew);
    nv++;    
  }

  // Enumerate all extreme points
  int origNumVars = originalInstance_->getNumVars();
  vector<int> S(origNumVars);
  for(int i = 0; i < origNumVars; i++) {
    S[i] = i;
  }
  vector<vector<double> > V;
  allExtreme(S, lb, ub, V);

  // Add lambda variables
  vector<VariablePtr> lambdavars;
  for(UInt i = 0; i < V.size(); i++) {
    VariablePtr vnew = VariablePtr(new Variable(nv,0.0,1.0,Continuous));
    variables_.push_back(vnew);
    lambdavars.push_back(vnew);
    nv++;    
  }    

  // Add the original linear constraints (on x)
  for(int i = 0; i < lcid.size(); i++) {
    const ConstraintPtr mlc = originalInstance_->getConstraint(lcid[i]);
    const LinearFunctionPtr olf = mlc->getLinearFunction();
    LinearFunctionPtr lf = LinearFunctionPtr(new LinearFunction());
    for(ConstVariableGroupIterator it = olf->varsBegin(); it != olf->varsEnd(); ++it) {
      lf->addTerm(xvars[it->first->getId()], it->second);
    }
    FunctionPtr f = (FunctionPtr) new Function(lf);
    ConstraintPtr c = (ConstraintPtr) new Constraint(f, mlc->getLb(), mlc->getUb());
    constraints_.push_back(c);    
#if defined(DEBUG_CONVEX_HULL)
      c->display();
#endif

  }

  // Write x in terms of extreme points
  for(int i = 0; i < origNumVars; i++) {
    LinearFunctionPtr lf = LinearFunctionPtr(new LinearFunction());
    lf->addTerm(xvars[i], -1.0);
    for(int j = 0; j < V.size(); j++) {
      lf->addTerm(lambdavars[j], V[j][i]);
    }
    //lf->display();
    FunctionPtr f = (FunctionPtr) new Function(lf);
    ConstraintPtr c = (ConstraintPtr) new Constraint(f, 0.0, 0.0);
    constraints_.push_back(c);    
    //XXX Should I set the ID?
#if defined(DEBUG_CONVEX_HULL)
      c->display();
#endif

  }

  // Write z in terms of extreme points
  for(int i = 0; i < mlcid.size(); i++) {
     LinearFunctionPtr lf = LinearFunctionPtr(new LinearFunction());
     lf->addTerm(zvars[i], -1.0);
     const ConstraintPtr mlc = originalInstance_->getConstraint(mlcid[i]);
     const FunctionPtr mlf = mlc->getFunction();

     for(int j = 0; j < V.size(); j++) {
       double zval = mlf->eval(V[j]);
       //cout << "zval = " << zval << endl;
       lf->addTerm(lambdavars[j], zval);
     }
     //lf->display();  cout << endl << endl;     
     FunctionPtr f = (FunctionPtr) new Function(lf);
     ConstraintPtr c = (ConstraintPtr) new Constraint(f, 0.0, 0.0);
     constraints_.push_back(c);    
#if defined(DEBUG_CONVEX_HULL)
      c->display();
#endif

  }

  // Add convexity constraint
  LinearFunctionPtr lf =  LinearFunctionPtr(new LinearFunction());
  for(int j = 0; j < V.size(); j++) {
    lf->addTerm(lambdavars[j], 1.0);
  }
  FunctionPtr f = (FunctionPtr) new Function(lf);
  ConstraintPtr c = (ConstraintPtr) new Constraint(f, 1.0, 1.0);
  constraints_.push_back(c);      

#if defined(DEBUG_CONVEX_HULL)
      c->display();
#endif


  LinearFunctionPtr olf = LinearFunctionPtr(new LinearFunction());

  // Add objective (on x vars only)
  const LinearFunctionPtr originalInstance_olf = originalInstance_->getObjective()->getLinearFunction();
  for (ConstVariableGroupIterator it = originalInstance_olf->varsBegin(); 
       it != originalInstance_olf->varsEnd(); ++it) {
    olf->addTerm(xvars[it->first->getId()], it->second);    
  }
  FunctionPtr of = (FunctionPtr) new Function(olf);
  objective_ = ObjectivePtr(new Objective(of, 0));    

}
Ejemplo n.º 4
0
void
MultilinearFormulation::makeRowByRow()
{

  // First we do some processing of the instance to determine how many 
  //  multilinear or quadratic constraints, and we store the indicies
  // for later    
  vector<int> lcid;
  vector<int> mlcid;
  for(ConstConstraintIterator it = originalInstance_->consBegin(); 
      it != originalInstance_->consEnd(); ++it) {
    FunctionType ft = (*it)->getFunctionType();
    if (ft == Multilinear) {
      mlcid.push_back((*it)->getId());
    }    
    else if (ft == Bilinear) {
      mlcid.push_back((*it)->getId());
    }
    else if (ft == Linear) {
      lcid.push_back((*it)->getId());
    }
  }
    
  // add x variables
  vector<double> lb;
  vector<double> ub;

  vector<VariablePtr> xvars;
  int nv = 0;
  for (ConstVariableIterator it = originalInstance_->varsBegin(); 
       it != originalInstance_->varsEnd(); ++it) {
    VariablePtr v = *it;
    VariablePtr vnew = VariablePtr(new Variable(nv, v->getLb(), v->getUb(), v->getType()));
    lb.push_back(v->getLb());
    ub.push_back(v->getUb());

    variableMapping_.insert(make_pair(vnew,v));    
    variables_.push_back(vnew);
    xvars.push_back(vnew);
    nv++;
  }

  // Add the linear constraints
  for(int i = 0; i < lcid.size(); i++) {
    const ConstraintPtr mlc = originalInstance_->getConstraint(lcid[i]);
    const LinearFunctionPtr olf = mlc->getLinearFunction();
    LinearFunctionPtr lf = LinearFunctionPtr(new LinearFunction());
    for(ConstVariableGroupIterator it = olf->varsBegin(); it != olf->varsEnd(); ++it) {
      lf->addTerm(xvars[it->first->getId()], it->second);
    }
    FunctionPtr f = (FunctionPtr) new Function(lf);
    ConstraintPtr c = (ConstraintPtr) new Constraint(f, mlc->getLb(), mlc->getUb());
    constraints_.push_back(c);    
  }

  vector<VariablePtr> zvars;
  vector<vector< VariablePtr> > lambdavars(mlcid.size());  
  for(UInt i = 0; i < mlcid.size(); i++) {
    
    // One 'z' var per row
    VariablePtr zvar = VariablePtr(new Variable(nv, -INFINITY, INFINITY, Continuous));
    variables_.push_back(zvar);
    zvars.push_back(zvar);
    nv++;
    
    // Determine what variables appear nonlinearly...
    const ConstraintPtr mlc = originalInstance_->getConstraint(mlcid[i]);
    vector<int> nlvars = nonlinearVarsInConstraint(mlc);

#if defined(DEBUG_ROW_BY_ROW)
    cout << "Vars in ml constraint " << i << " have ix: ";
    for(vector<int>::iterator it = nlvars.begin(); it != nlvars.end(); ++it) {
      cout << (*it) << " ";
    }
    cout << endl;
#endif

    
    // Add lambda vars for this row
    vector<vector<double> > V;
    allExtreme(nlvars, lb, ub, V);
    for(UInt j = 0; j < V.size(); j++) {
      VariablePtr vnew = VariablePtr(new Variable(nv,0.0,1.0,Continuous));
      variables_.push_back(vnew);
      lambdavars[i].push_back(vnew);
      nv++;
    }

    for(UInt k = 0; k < nlvars.size(); k++) {
      // Create 'x' portion of this constraint
      LinearFunctionPtr lf = LinearFunctionPtr(new LinearFunction());
      lf->addTerm(xvars[nlvars[k]], -1.0);

      for(UInt j = 0; j < V.size(); j++) {
        lf->addTerm(lambdavars[i][j], V[j][k]);
      }                
      
      FunctionPtr f = (FunctionPtr) new Function(lf);
      ConstraintPtr c = (ConstraintPtr) new Constraint(f, 0.0, 0.0);
      constraints_.push_back(c);  

#if defined(DEBUG_ROW_BY_ROW)
      c->display();
#endif
    }
          
    
    // Now add the 'z' var definition
    LinearFunctionPtr tlf = LinearFunctionPtr(new LinearFunction());
    tlf->addTerm(zvar, -1.0);
    const FunctionPtr mlf = mlc->getFunction();
    const LinearFunctionPtr mlf_lf = mlf->getLinearFunction();
    

    for(int j = 0; j < V.size(); j++) {
      //XXX Kludgy, but just create a vector for the (full) point
      int onv = originalInstance_->getNumVars();
      vector<double> xe(onv,0.0);
      for(UInt k = 0; k < nlvars.size(); k++) {
        xe[nlvars[k]] = V[j][k];
      }
      double zval = mlf->eval(xe);
      // Need to subtract off linear part, since yuo keep those variables in
      zval -= mlf_lf->eval(xe);

      tlf->addTerm(lambdavars[i][j], zval);
    }
    FunctionPtr tf = (FunctionPtr) new Function(tlf);
    ConstraintPtr tc = (ConstraintPtr) new Constraint(tf, 0.0, 0.0);
    constraints_.push_back(tc);  

#if defined(DEBUG_ROW_BY_ROW)
      tc->display();
#endif
    
    // Now add the linear constraint involving linear x and t
    LinearFunctionPtr xtlf = LinearFunctionPtr(new LinearFunction());
    for(ConstVariableGroupIterator it = mlf_lf->varsBegin(); it != mlf_lf->varsEnd(); ++it) {
      xtlf->addTerm(xvars[it->first->getId()], it->second);
    }
    xtlf->addTerm(zvar,1.0);
    FunctionPtr xtf = (FunctionPtr) new Function(xtlf);
    ConstraintPtr xtc = (ConstraintPtr) new Constraint(xtf, mlc->getLb(), mlc->getUb());
    constraints_.push_back(xtc);  

#if defined(DEBUG_ROW_BY_ROW)
    xtc->display();
#endif


    // Also add sum (lambda) = 1
    LinearFunctionPtr convex_lf =  LinearFunctionPtr(new LinearFunction());
    for(int j = 0; j < V.size(); j++) {
      convex_lf->addTerm(lambdavars[i][j], 1.0);
    }
    FunctionPtr convex_f = (FunctionPtr) new Function(convex_lf);
    ConstraintPtr convex_c = (ConstraintPtr) new Constraint(convex_f, 1.0, 1.0);
    constraints_.push_back(convex_c);      

#if defined(DEBUG_ROW_BY_ROW)
    convex_c->display();
#endif

    
  }



  LinearFunctionPtr olf = LinearFunctionPtr(new LinearFunction());
  const LinearFunctionPtr originalInstance_olf = originalInstance_->getObjective()->getLinearFunction();
  for (ConstVariableGroupIterator it = originalInstance_olf->varsBegin(); 
       it != originalInstance_olf->varsEnd(); ++it) {
    olf->addTerm(xvars[it->first->getId()], it->second);    
  }
  FunctionPtr of = (FunctionPtr) new Function(olf);
  objective_ = ObjectivePtr(new Objective(of, 0));    
  

}