Multiple Tableau Gomory Cuts


Detailed Description

This code present an implementation of some of the ideas on generating cuts from multiple tableau rows at once. It has a basic cut-separation code, and a cut-callback implementation for CPLEX


Files

file  mt_cplex.ex.c
file  mt_cplex_cbk.c
file  mt_cplex_cbk.h
file  mt_cutpool.c
file  mt_cutpool.h
file  mt_gomory.c
file  mt_gomory.h
file  mt_SL.c
file  mt_SL.h
file  mt_T1.c
file  mt_T1.h
file  mt_T2.c
file  mt_T2.h
file  mt_tableau.c
file  mt_tableau.h
file  mt_types.h

Data Structures

struct  MTccbk_info_t
 structure holding basic calllback information, used to control callback behavior and statistic More...
struct  MTcut_t
 structure to hold a cut More...
struct  MTcutHeap_t
 simple structure (based on heaps) that holds up to certain number of cuts, and keep (acording to a criteria) the best ones. More...
struct  MTcutHeapCn_t
 basic conector for the cut heap More...
struct  MTcutHeapControl_t
 configuration strucutre for MTcutHeap_t More...
struct  MTgomory_ccbk_t
 structure holding the working information for the all callbacks function More...
struct  MTlp_t
 formulation description for a problem More...
struct  MTns_ccbk_t
 information for node solve callback More...
struct  MTrowmatrix_t
 structure holding a compresed matrix in row form More...
struct  MTsol_t
 structure to hold a solution to an LP More...

#MTccbk_info_t set/get functions

#define MTcbkiSetCuts(__cbki, __val)   ((__cbki)->use_cuts=(__val))
#define MTcbkiSetCutsAtRoot(__cbki, __val)   ((__cbki)->root_only=(__val))
#define MTcbkiSetCutStyle(__cbki, __val)   ((__cbki)->cut_style=(__val))
#define MTcbkiSetMaxNodeCuts(__cbki, __val)   ((__cbki)->max_cuts=(size_t)(__val))
#define MTcbkiSetMaxRoundCuts(__cbki, __val)   ((__cbki)->max_cpr=(__val))
#define MTcbkiSetMaxRows(__cbki, __val)   ((__cbki)->max_rows=(__val))

Profiling information for the cuts

#define MTccbk_display_sumary(__fout)
 display summary information for the call-back calls
double MTccbk_discarded_ratio
 minimum ratio of discarded inequality
double MTccbk_discarded_violation
 maximum violation of discarded inequality
int MTccbk_fail_ratio
 Number of discarded cuts because high ratio.
int MTccbk_fail_tableau_btype
 Type of basic variable for tableau is not integer.
int MTccbk_fail_tableau_coeff
 Number of discarded cuts because tableaus has wrong coefficient for basic variable.
int MTccbk_fail_tableau_frac
 Number of discarded cuts because fraction is too small.
int MTccbk_fail_tableau_nbasic
 Number of discarded cuts because tableaus has wrong number of basic variables.
int MTccbk_fail_tableau_ratio
 display summary information for the call-back calls
int MTccbk_fail_violation
 Number of discarded cuts because low violation.

#MTgomory_ccbk_t set/get functions

#define MTccbkSetCutDominance(__cbdata, __val)   MTcutHeapSetDominance((&((__cbdata)->cuth)),__val)
#define MTccbkSetCuts(__cbdata, __val)   MTcbkiSetCuts((&((__cbdata)->info)),__val)
#define MTccbkSetCutsAtRoot(__cbdata, __val)   MTcbkiSetCutsAtRoot((&((__cbdata)->info)),__val)
#define MTccbkSetCutSel(__cbdata, __val)   MTcutHeapSetCS((&((__cbdata)->cuth)),__val)
#define MTccbkSetCutStyle(__cbdata, __val)   MTcbkiSetCutStyle((&((__cbdata)->info)),__val)
#define MTccbkSetMaxNodeCuts(__cbdata, __val)   MTcbkiSetMaxNodeCuts((&((__cbdata)->info)),__val)
#define MTccbkSetMaxRoundCuts(__cbdata, __val)
#define MTccbkSetMaxRows(__cbdata, __val)   MTcbkiSetMaxRows((&((__cbdata)->info)),__val)

static parameters for the main program

static int cut_dominance = 0
 cut dominance check enable
static MTcutSelection_t cut_selection = MTcsDefault
 cut selection strategy
static int cut_style = 1
 cut generation strategy
static int cuts_at_root_only = 1
 cuts at root only
static CPXENVptr env = 0
 global cplex environment
static CPXLPptr lp = 0
 global lp pointer
static int max_cuts_per_node = 1000
 max cuts per node
static int max_cuts_per_round = 10
 max cuts per round
static double max_rtime = INT_MAX
 maximum running time
static int max_tableau_rows = 5
 maximum number of tableau-rows to use
static unsigned long memlimit = UINT_MAX
 maximum memory usage
static char options [1024] = "cplex.opt"
 CPLEX option file name.
static char problem [1024] = ""
 input problem file
static int use_cuts = 1
 use callback cuts
static int use_solve_cbk = 0
 use solve callback

Defines

#define __MT_DOM
 helper to set dominance
#define __MT_NS_LVL   0
#define MT_CCBK_MAX_RATIO   0x1p15
 Maximum allowed cut ratio. The cut ratio is defined as the ratio between the maximum absolute value and minimum absolute value for the non-zero coefficients. It is importatn because is a measure of the error inducing capabilities of the inequality.
#define MT_CCBK_MIN_FRAC   0x1p-12
 Minimum fractionality required for a tableau to be considered, i.e. if abs(x-round(x))< min_frac we discard the asociated tableau.
#define MT_CCBK_MIN_TABLEAU_RATIO   0x1p-12
 Minimum tableau row ratio. This is defined as minabsval/maxabsval for a tableau row. Tableau rows that have a ratio bellow the given threshold will not be considered at the moment of cut genweration.
#define MT_CCBK_MIN_VIO   0x1p-10
 Minimum violation for a cut to be added. Note that we scale the cut in such a way that the maximum absolute value of a non-zero coefficient is one.
#define MT_CCBK_USE_NAMES   1
 if se to one, use real variable names while displaying cut/tableau information
#define MT_CPLEX_SENSE   'G'
#define MT_DEBUG_LVL   100
 debug level for the module
#define MT_VERB_LVL   100
 verbose level for the module
#define MT_ZERO_TOL   0x1p-20
 tolerance level for zero, i.e. all absolute values bellow MT_ZERO_TOL are considered as zero.
#define MTccbk_info_init(__cbinf)   memset(__cbinf,0,sizeof(MTccbk_info_t))
 initialize an MTccbk_info_t structure
#define MTcplexCHECKRVALG(__env, __rval, __where)
 internal debugging information for CPLEX
#define MTcut_init(__cut)   memset(__cut,0,sizeof(MTcut_t))
 initialize a cut structure
#define MTcutHeap_get_free_cut(__cuth)   ((__cuth)->next_free)
 return a pointer to the next free MTcut_t structure in the set
#define MTcutHeap_get_ncuts(__cuth)   (__cuth->heap.sz)
 return the number of stored cuts in the given MTcutHeap_t structure.
#define MTcutHeapCnClear(__chcn)
 clear all internally allocated memory
#define MTcutHeapCnInit(__chcn)
 initialize an MTcutHeapCn_t structure
#define MTcutHeapSetCS(__cuth, __cutsel)   ((__cuth)->control.cut_score = (__cutsel))
 set the cutselection strategy
#define MTcutHeapSetDominance(__cuth, __dom)   ((__cuth)->control.dominance_check = (__dom))
 set the dominance strategy
#define MTgomory_ccbk_clear(__cbdata)
 free any internally allocated memory inside an MTgomory_ccbk_t structure.
#define MTgomory_ccbk_init(__cbdata)
 initialize an MTgomory_ccbk_t structure.
#define MTgomory_ccbk_rsz(__cbdata, __nrowsz, __ncolsz, __ntbsz)
 if required (i.e. __nrowsz > nrowsz || __ncolsz > colsz || __ntbsz > tbsz) , resize an MTgomory_ccbk_t structure
#define MTlp_init(__lp)   memset((__lp),0,sizeof(MTlp_t))
 initialize an MTlp_t structure
#define MTns_ccbk_clear(__ns_ccbk)
 clear any internally allocated memory within an MTns_ccbk_t structure.
#define MTns_ccbk_display(__ns_ccbk)
 Display stats in the given MTns_ccbk_t structure.
#define MTns_ccbk_init(__ns_ccbk)
 initialize an MTns_ccbk_t structure
#define MTns_ccbk_rsz(__ns_ccbk, __ncsz, __nrsz)
 resize node-solve arrays
#define MTrowmatrix_init(__mt)   memset((__mt),0,sizeof(MTrowmatrix_t))
 Initialize an MTrowmatrix_t structure too an empty matrix.
#define MTsol_init(__sol)   memset(__sol,0,sizeof(MTsol_t))
 initialize a MTsol_t structure
#define UPDATE_F_AND_VAL(__f, __cstat, __bound, __val)
 simple handler to update f according to the complementation of variables

Enumerations

enum  MTcutSelection_t {
  MTcsRAW, MTcsVIOL, MTcsINVLINF, MTcsINVL2,
  MTcsDefault = MTcsRAW
}
 types of cut score functions More...

Functions

int main (int argc, char **argv)
 this function reads a problem file, set the cut-callback, and call the MIP solver
static int mem_limits (void)
 function to handle resource usage limits
void MTccbk_info_display (const MTccbk_info_t *const info, FILE *out)
 display callback info
int MTccbk_info_process (MTccbk_info_t *const info, CPXCENVptr env, void *cbdata, int wherefrom, int *const action)
 test if we should proceed with the cuting process, and update information
int MTcheckTabrow (const MTlp_t *const lp, const double *tableau, int *const status, double *const ratio, int *const etnz)
 given a tableau in extended notation, verify that it is indeed a numerically stable and correct tableau.
int MTcompressTabRow (const MTlp_t *const lp, const MTsol_t *const sol, double *rowval, int *rowind, double *const rhs, int *const nz)
 given a tableau row in extended format, store it as a compresed form, it include setting row to row - integer_part[row] and updating the rhs.
void MTcompute_f (const MTlp_t *const lp, const double *const tableau, const double *const x, const double *const slack, double *const f)
 compute ax = f for a tableau
int MTcompute_fractional_vars (const MTlp_t *const lp, const double *const x, int *const fracvars, int *const nfrac, const double minitgap)
 get the indices of fractional (structural) variables
int MTcompute_integer_vars (const MTlp_t *const lp, const double *const x, int *const intvars, int *const nint, const double intgap)
 get the indices of (almost) integer variables
int MTcplex_binv_to_tableau (const MTlp_t *const lp, const double *const binvrow, double *const tableau)
 given an LP desvription, and CPLEX binvrow, compute the actual tableau (including logicals, which are in position i+ncols of the tableau for the i-th logical/slack variable, but we force that the coefficient of equality logicals is zero)
void MTcplex_display_compress_tableau (FILE *file, const MTlp_t *const lp, const double *const matval, const int *const matind, const int nz)
 display a compress tableau into the given file
void MTcplex_display_tableau (FILE *file, const MTlp_t *const lp, const double *const tableau)
 display a tableau into the given file
int MTcplex_test_cut (const MTcut_t *const cut, const MTlp_t *const lp)
 test if a given cut is valid to the given cplex LP
int MTcut_check_dominance (const MTcut_t *const cut1, const MTcut_t *const cut2)
 check for dominance between cuts.
void MTcut_clear (MTcut_t *const MTcut)
 free any internal memory asociated with an MTcut_t structure
void MTcut_copy (MTcut_t *const dest, const MTcut_t *const src)
 copy a cut from one structure to the next
void MTcut_rsz (MTcut_t *const MTcut, const int MTnsz)
 resize a cut structure
int MTcutHeap_add_cut (const MTlp_t *const lp, MTcutHeap_t *const cuth)
 add the cut (obtained by MTcutHeap_get_free_cut) to the structure, and possibly discard a worst cut from the heap
MTcut_tMTcutHeap_pop_cut (MTcutHeap_t *const cuth)
 pop the worst remaining cut within an MTcutHeap_t structure. Note that the memory asociated with the returned MTcut_t structure is handled by the MTcutHeap_t structure, and should not be freed by the user.
void MTcutHeapClear (MTcutHeap_t *const cuth)
 Clear any internal memory asociated with an MTcutHeap_t structure. Note that after a call to this function a new call to MTcutHeapInit should be made before using any other function over the strucutre.
static int MTcutHeapCn_compute_val (MTcutHeap_t *const cuth, MTcutHeapCn_t *const chcn)
 given a MTcutHeapCn_t compute the value asociated with the stored cut and store it in MTcutHeapCn_t::hcn, inside dbl_EGeHeapCn_t::val, the function used to evaluate depend on the configuration for the cut heap. The smallest the value, the worst the cut is.
void MTcutHeapInit (MTcutHeap_t *const cuth, size_t sz)
 initialize an MTcutHeap_t structure
void MTcutNzAlloc (void)
void MTcutNzFree (void)
void MTdisplay_lp (const MTlp_t *const lp, const MTsol_t *const sol, FILE *stream)
 display an LP to the given file
int MTget_best_k_integer_tableau_rows (const MTlp_t *const lp, const MTsol_t *const sol, const int max_tableau, CPXCENVptr env, CPXLPptr CPXlp, MTrowmatrix_t *const tb, double *const f)
 Get the best K tableau (asociated to integer basic variables), note that we check for wrongly computed tableau. Note that we compute Z = f + 0.5+Ax, where A is the tableau and f the fractional value (displaced by 1/2), moreover, all basic variables have zero coefficient in the resulting tableau, and all variables are complemented to their bounds.
int MTget_best_k_tableau_rows (const MTlp_t *const lp, const MTsol_t *const sol, const int max_tableau, CPXCENVptr env, CPXLPptr CPXlp, MTrowmatrix_t *const tb, double *const f)
 Get the best K tableau (asociated to fractional variables), note that we check for wrongly computed tableau. Note that we compute Z = f + 0.5+Ax, where A is the tableau and f the fractional value (displaced by 1/2), moreover, all basic variables have zero coefficient in the resulting tableau, and all variables are complemented to their bounds.
int MTget_sol (const MTlp_t *const lp, CPXCENVptr env, CPXLPptr CPXlp, void *cbdata, int wherefrom, MTsol_t *const sol)
 obtain a solution to the current node LP
int MTgomory_ccbk (CPXCENVptr env, void *cbdata, int wherefrom, void *cbhandle, int *useraction_p)
 Interface for the CPLEX-callback for the MTgomoryCut function.
int MTgomoryCut (const int nrows, const int ncols, const double *const rowval, const int *const rowind, const int *const rowbeg, const double *const f, const double *const base, double *const cutval, double *const work)
 Given a set of tableau rows in column format, return the asociated multi-tableau gomory cut. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols. The returned constraint is $\displaystyle \sum\limits_{i\in I}s_ib_i\geq 1$ where $ \displaystyle b_i=$ cutval[i], $\displaystyle \forall i\in I$. And where the set used to cut is $\displaystyle L=\left\{x:c\sigma x\leq c_o ,\forall\sigma\in\mathrm{diag}\left(\{-1,1\}^{m}\right)\right\}$ and $\displaystyle c=$ base, $\displaystyle c_o = \frac12\sum\limits_{i=1}^mc_i$, where we assume that $c_i\neq 0\forall i=1,\ldots,m$. Finally, note that $\displaystyle b_i=\max\limits_{\sigma:c\sigma r_i>0}\left\{\frac{c\sigma r_i}{c_o-c\sigma f}\right\}$.

void MTintPermSort (const int sz, int *const perm, const int *const elem)
 permutation sort for integers
void MTlp_clear (MTlp_t *const MTlp)
 free any internally allocated memory inside an MTlp_t structure.
int MTlp_load_problem (CPXCENVptr env, MTlp_t *const lp, CPXLPptr CPXlp, void *cbdata, int wherefrom)
 Given an LP problem pointer, load it into an (initialized) MTlp_t structure.
static void MTlp_rsz (MTlp_t *const MTlp, const int rsz, const int csz, const int tsz)
 resize an lp structure
int MTnode_solve (CPXCENVptr env, void *cbdata, int wherefrom, void *cbhandle, int *useraction_p)
 simple node solve callback, we use it to count primal degenerancy, i.e. the number of variables with zero reduced cost, this should indicate posibility of doing nested objective funcions to avoid primal degenerancy
int MTraw_to_complemented_cut (const MTlp_t *const lp, const MTsol_t *const sol, double *const extended, MTcut_t *const cut)
 transform an extened non-complemented cut to a compresed and complemented form. We assume the extended cut is of the form ax >= 1. Note that we expect no coefficient for basic variables in the cut, moreover, the rhs of the inequality is assumed to be 1 (in raw format).
int MTraw_to_uncomplemented_cut (const MTlp_t *const lp, const double *const extended, const int nactive, const int *const active, MTcut_t *const cut, int *const status)
 transform an extened non-complemented cut to a compresed for. We assume the extended cut is of the form ax >= 1.m Note that we expect no coefficient for basic variables in the cut.
int MTread_cplex_options (CPXENVptr env, FILE *inputfile)
 read optionfile and set parameter to CPLEX
int MTreadRTableau (FILE *input, int *const nrows, int *const ncols, double **const rowval, int **const rowind, int **const rowbeg, double **const f)
 read a tableau in row format from a file (see MTgomoryCut for details on the tableau information).
void MTrowmatrix_clear (MTrowmatrix_t *const MTrm)
 free any internally allocated memory inside an MTrowmatrix_t structure.
void MTrowmatrix_rsz (MTrowmatrix_t *const MTrm, const int MTrsz, const int MTnsz)
 resize a MTrowmatrix_t structure
int MTsl1 (const double *const ineq, MTrowmatrix_t *const tb, double *const work, const char *const vtype, const int mode)
 compute in work the value ineq * r_i as described by the given inequality and tableau, we assume that work is zero, moreover, given the variable status (interger/binary) and an operation mode, try to do nothing (mode 0), try to get ineq*r[i] >0 (mode 1) or try to get ineq*r[i] < 0 (mode 2), note that mode != 0 will modify the tableau
int MTsl2 (const double *const ineq, double *const cutval, double *const work, const int *const active, const int nactive)
 , set cutval[i] = max(cutval[i],work[i] > 0 ? work[i]/ineq[3]: 0); and reset work to zero
int MTsl3 (const double *const f, double *const ineq, const int dim)
 given an inequality and f compute ineq[dim+1] and (possible change signs) to ensure that ineq[dim+1] > 0;
int MTsl4 (double *const ineq, const double x0, const double _y0, const double x1, const double _y1)
 compute the line that pases through the two given points
int MTsl5 (int *const nactive, int *const active, const int nvars, const int nz, const double *const compval, const int *const compind)
 compute active set from the given compresed vector set
int MTsl_ccbk (CPXCENVptr env, void *cbdata, int wherefrom, void *cbhandle, int *useraction_p)
 Interface for the CPLEX-callback for the MTt1Cut function.
int MTslCut (MTgomory_ccbk_t *const data, CPXCENVptr env, CPXLPptr lp)
 The idea of these cuts is to use some ideas from the paper ``On the relative strengt of split, triangle and cuadrilateral cuts'' from Cornuejols et. al. The main point is that splits are bad (compared against triangles or quadrilaterals) whenever we have (an almost) integer part in $f_i$ for some coordiante $i$, and a wide base, this can be seen on the following figure: $ \hspace*{\textwidth} \vspace*{2in} x_1 = 1 \psset{arrows=->,fillstyle=solid,fillcolor=black,linecolor=black} \begin{pspicture}(0,0)(5,3) \multido{\i=0+1}{6}{{ \pscircle(\i,0){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,1){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,2){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,3){1mm} }} \rput(2.5,2){$f$} \psline(2.5,2)(0,.5)\nbput{$r_1$} \psline(2.5,2)(4,1)\nbput{$r_2$} \psline[arrows=|-|](1,1)(4,1)\nbput{M} \end{pspicture} $.
void MTsol_clear (MTsol_t *const MTsol)
 clear any internal memory asociated with an MTsol_t structure
void MTsol_rsz (MTsol_t *const MTsol, const int cnz, const int rnz)
 resize a solution structure
int MTt1_ccbk (CPXCENVptr env, void *cbdata, int wherefrom, void *cbhandle, int *useraction_p)
 Interface for the CPLEX-callback for the MTt1Cut function.
int MTt1Cut (const MTrowmatrix_t *const tb, MTlp_t *const lp, MTcutHeap_t *const cuth, const double *const f, double *const extcut, double *const work, const int cut_style)
 compute the best cuts found using T1 sets as a base, where

\[ T1_n=\left\{ x\in\mathbb{R}^n: x_i\geq 0, \sum\limits_{i=1}^n x_i\leq n \right\},\]

and where we try all possible $2^n$ reorientation of the axis. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols.

int MTt2_ccbk (CPXCENVptr env, void *cbdata, int wherefrom, void *cbhandle, int *useraction_p)
 Interface for the CPLEX-callback for the MTt2Cut function.
int MTt2Cut (const MTrowmatrix_t *const tb, MTlp_t *const lp, MTcutHeap_t *const cuth, const double *const f, double *const extcut, double *const work, const int cut_style)
 compute the best cuts found using T2 sets as a base, where

\[ T2_n=\left\{x\in\mathbb{R}^n: \begin{matrix} R1:&-x_1 & \leq 0\\ R2:& x_1-x_2 & \leq 1\\ R3:& x_1+x_2-x_3 & \leq 2\\ \vdots & \vdots & \vdots\\ Rn:& x_1+\ldots+x_{n-1}-x_n&\leq n-1\\ R(n+1):& x_1+\ldots+x_{n-1}+x_n&\leq n\\ \end{matrix}\right\} \]

, and where we try all possible $n!2^n$ reorientation of the axis and permutation of the rows. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols.

int MTuncomplemented_to_complemented_cut (const MTlp_t *const lp, const MTsol_t *const sol, MTcut_t *const cut)
 transform an uncomplemented cut to a complemented form. Note that we expect no coefficient for basic variables in the cut.
int MTwriteRTableau (FILE *out, const int nrows, const int ncols, const double *const rowval, const int *const rowind, const int *const rowbeg, const double *const f)
 write a tableau in row format into a file (see MTgomoryCut for details on the tableau information).
int parseargs (int argc, char **argv)
 parsing parameters function, also read problem and options file
static void sighandler (int s)
void usage (char *const prog)
 usage function

Variables

char __mt_cplex_errbuf [4096]
 internal error string for cplex
char __mt_cplex_errbuf [4096]
 internal error string for cplex
double MTccbk_discarded_ratio = DBL_MAX
 minimum ratio of discarded inequality
double MTccbk_discarded_violation = 0
 maximum violation of discarded inequality
int MTccbk_fail_ratio = 0
 Number of discarded cuts because high ratio.
int MTccbk_fail_tableau_btype = 0
 Type of basic variable for tableau is not integer.
int MTccbk_fail_tableau_coeff = 0
 Number of discarded cuts because tableaus has wrong coefficient for basic variable.
int MTccbk_fail_tableau_frac = 0
 Number of discarded cuts because fraction is too small.
int MTccbk_fail_tableau_nbasic = 0
 Number of discarded cuts because tableaus has wrong number of basic variables.
int MTccbk_fail_tableau_ratio = 0
int MTccbk_fail_violation = 0
 Number of discarded cuts because low violation.
static int * MTcut1_nz = 0
 internal memory for dominance check
static int * MTcut2_nz = 0
const char MTversion [1024]
 Version String of the program.


Define Documentation

#define __MT_DOM
 

Value:

do{\
  if(v2 + MT_ZERO_TOL < v1) \
  {\
    if(domset && dom != -1) return 0;\
    domset=1;\
    dom=-1;\
    equal = 0; \
  }\
  else if( v1 + MT_ZERO_TOL < v2)\
  {\
    if(domset && dom != 1) return 0;\
    domset=1;\
    dom=1;\
    equal = 0;\
  }\
  }while(0)
helper to set dominance

Definition at line 176 of file mt_cutpool.c.

#define __MT_NS_LVL   0
 

Definition at line 1847 of file mt_cplex_cbk.c.

#define MT_CCBK_MAX_RATIO   0x1p15
 

Maximum allowed cut ratio. The cut ratio is defined as the ratio between the maximum absolute value and minimum absolute value for the non-zero coefficients. It is importatn because is a measure of the error inducing capabilities of the inequality.

Definition at line 48 of file mt_cplex_cbk.h.

#define MT_CCBK_MIN_FRAC   0x1p-12
 

Minimum fractionality required for a tableau to be considered, i.e. if abs(x-round(x))< min_frac we discard the asociated tableau.

Definition at line 42 of file mt_cplex_cbk.h.

#define MT_CCBK_MIN_TABLEAU_RATIO   0x1p-12
 

Minimum tableau row ratio. This is defined as minabsval/maxabsval for a tableau row. Tableau rows that have a ratio bellow the given threshold will not be considered at the moment of cut genweration.

Definition at line 53 of file mt_cplex_cbk.h.

#define MT_CCBK_MIN_VIO   0x1p-10
 

Minimum violation for a cut to be added. Note that we scale the cut in such a way that the maximum absolute value of a non-zero coefficient is one.

Definition at line 38 of file mt_cplex_cbk.h.

#define MT_CCBK_USE_NAMES   1
 

if se to one, use real variable names while displaying cut/tableau information

Definition at line 61 of file mt_cplex_cbk.h.

#define MT_CPLEX_SENSE   'G'
 

if we have CPLEX 11 or later, we use sense as 'g' instead of 'G', because, in theory, it enables cutpool management, we will see

Definition at line 542 of file mt_cplex_cbk.h.

#define MT_DEBUG_LVL   100
 

debug level for the module

Definition at line 99 of file mt_gomory.h.

#define MT_VERB_LVL   100
 

verbose level for the module

Definition at line 95 of file mt_gomory.h.

#define MT_ZERO_TOL   0x1p-20
 

tolerance level for zero, i.e. all absolute values bellow MT_ZERO_TOL are considered as zero.

Definition at line 104 of file mt_gomory.h.

#define MTcbkiSetCuts __cbki,
__val   )     ((__cbki)->use_cuts=(__val))
 

Definition at line 168 of file mt_cplex_cbk.h.

#define MTcbkiSetCutsAtRoot __cbki,
__val   )     ((__cbki)->root_only=(__val))
 

Definition at line 167 of file mt_cplex_cbk.h.

#define MTcbkiSetCutStyle __cbki,
__val   )     ((__cbki)->cut_style=(__val))
 

Definition at line 172 of file mt_cplex_cbk.h.

#define MTcbkiSetMaxNodeCuts __cbki,
__val   )     ((__cbki)->max_cuts=(size_t)(__val))
 

Definition at line 170 of file mt_cplex_cbk.h.

#define MTcbkiSetMaxRoundCuts __cbki,
__val   )     ((__cbki)->max_cpr=(__val))
 

Definition at line 171 of file mt_cplex_cbk.h.

#define MTcbkiSetMaxRows __cbki,
__val   )     ((__cbki)->max_rows=(__val))
 

Definition at line 169 of file mt_cplex_cbk.h.

#define MTccbk_display_sumary __fout   ) 
 

Value:

do{\
  fprintf(__fout,"MTccbk sumary:\n");\
  fprintf(__fout,"\tDiscarded tableau rows: %4d (coeff %d, nbasic %d,"\
          " ratio %d, type %d, fraction %d)\n", MTccbk_fail_tableau_nbasic + \
          MTccbk_fail_tableau_coeff + MTccbk_fail_tableau_ratio + \
          MTccbk_fail_tableau_btype + MTccbk_fail_tableau_frac,\
          MTccbk_fail_tableau_coeff, MTccbk_fail_tableau_nbasic, \
          MTccbk_fail_tableau_ratio, MTccbk_fail_tableau_btype, \
          MTccbk_fail_tableau_frac);\
  fprintf(__fout,"\tDiscarded by ratio:     %4d (next %.5lg)\n",\
          MTccbk_fail_ratio, MTccbk_discarded_ratio);\
  fprintf(__fout,"\tDiscarded by violation: %4d (next %.5lg)\n",\
          MTccbk_fail_violation, MTccbk_discarded_violation);\
  }while(0)
display summary information for the call-back calls

Parameters:
__fout where to write the sumary

Definition at line 224 of file mt_cplex_cbk.h.

#define MTccbk_info_init __cbinf   )     memset(__cbinf,0,sizeof(MTccbk_info_t))
 

initialize an MTccbk_info_t structure

Definition at line 176 of file mt_cplex_cbk.h.

#define MTccbkSetCutDominance __cbdata,
__val   )     MTcutHeapSetDominance((&((__cbdata)->cuth)),__val)
 

Definition at line 242 of file mt_cplex_cbk.h.

#define MTccbkSetCuts __cbdata,
__val   )     MTcbkiSetCuts((&((__cbdata)->info)),__val)
 

Definition at line 245 of file mt_cplex_cbk.h.

#define MTccbkSetCutsAtRoot __cbdata,
__val   )     MTcbkiSetCutsAtRoot((&((__cbdata)->info)),__val)
 

Definition at line 244 of file mt_cplex_cbk.h.

#define MTccbkSetCutSel __cbdata,
__val   )     MTcutHeapSetCS((&((__cbdata)->cuth)),__val)
 

Definition at line 243 of file mt_cplex_cbk.h.

#define MTccbkSetCutStyle __cbdata,
__val   )     MTcbkiSetCutStyle((&((__cbdata)->info)),__val)
 

Definition at line 246 of file mt_cplex_cbk.h.

#define MTccbkSetMaxNodeCuts __cbdata,
__val   )     MTcbkiSetMaxNodeCuts((&((__cbdata)->info)),__val)
 

Definition at line 248 of file mt_cplex_cbk.h.

#define MTccbkSetMaxRoundCuts __cbdata,
__val   ) 
 

Value:

do{\
  MTgomory_ccbk_t*const __MTcbdata = (__cbdata);\
  MTcbkiSetMaxRoundCuts((&(__MTcbdata->info)),__val);\
  MTcutHeapClear(&(__MTcbdata->cuth));\
  MTcutHeapInit(&(__MTcbdata->cuth),((size_t)(__MTcbdata->info.max_cpr)));}while(0)

Definition at line 249 of file mt_cplex_cbk.h.

#define MTccbkSetMaxRows __cbdata,
__val   )     MTcbkiSetMaxRows((&((__cbdata)->info)),__val)
 

Definition at line 247 of file mt_cplex_cbk.h.

#define MTcplexCHECKRVALG __env,
__rval,
__where   ) 
 

Value:

do{\
  if(__rval)\
  {\
    CPXgeterrorstring(__env, __rval, __mt_cplex_errbuf);\
    fprintf(stderr,__mt_cplex_errbuf);\
    __EG_PRINTLOC__;\
    goto __where;\
  }\
}while(0)
internal debugging information for CPLEX

Parameters:
__rval cplex return value (zero means success)
env as returned by CPXopenCPLEX function.
where to jump if an error is found.

Definition at line 68 of file mt_cplex_cbk.h.

#define MTcut_init __cut   )     memset(__cut,0,sizeof(MTcut_t))
 

initialize a cut structure

Parameters:
__cut pointer to an MTcut_t structure

Definition at line 38 of file mt_cutpool.h.

#define MTcutHeap_get_free_cut __cuth   )     ((__cuth)->next_free)
 

return a pointer to the next free MTcut_t structure in the set

Parameters:
__cuth structure to initialize

Definition at line 71 of file mt_cutpool.h.

#define MTcutHeap_get_ncuts __cuth   )     (__cuth->heap.sz)
 

return the number of stored cuts in the given MTcutHeap_t structure.

Parameters:
__cuth structure where we perform the operation
Returns:
the number of cuts in the given MTcutHeap_t structure.

Definition at line 107 of file mt_cutpool.h.

#define MTcutHeapCnClear __chcn   ) 
 

Value:

do{\
  MTcutHeapCn_t*const __MTchcn = (__chcn);\
  dbl_EGeHeapCnClear(&(__MTchcn->hcn));\
  MTcut_clear(&(__MTchcn->cut));}while(0)
clear all internally allocated memory

Definition at line 57 of file mt_cutpool.h.

#define MTcutHeapCnInit __chcn   ) 
 

Value:

do{\
  MTcutHeapCn_t*const __MTchcn = (__chcn);\
  dbl_EGeHeapCnInit(&(__MTchcn->hcn));\
  MTcut_init(&(__MTchcn->cut));}while(0)
initialize an MTcutHeapCn_t structure

Definition at line 51 of file mt_cutpool.h.

#define MTcutHeapSetCS __cuth,
__cutsel   )     ((__cuth)->control.cut_score = (__cutsel))
 

set the cutselection strategy

Definition at line 66 of file mt_cutpool.h.

#define MTcutHeapSetDominance __cuth,
__dom   )     ((__cuth)->control.dominance_check = (__dom))
 

set the dominance strategy

Definition at line 63 of file mt_cutpool.h.

#define MTgomory_ccbk_clear __cbdata   ) 
 

Value:

do{\
  MTgomory_ccbk_t*const __ccbk = (__cbdata);\
  EGfree(__ccbk->work);\
  EGfree(__ccbk->extcut);\
  EGfree(__ccbk->base);\
  EGfree(__ccbk->f);\
  __ccbk->rowsz = __ccbk->colsz = 0;\
  MTsol_clear(&(__ccbk->sol));\
  MTcut_clear(&(__ccbk->cut));\
  MTrowmatrix_clear(&(__ccbk->tb));\
  MTcutHeapClear(&(__ccbk->cuth));\
  MTlp_clear(&(__ccbk->lp));}while(0)
free any internally allocated memory inside an MTgomory_ccbk_t structure.

Parameters:
__cbdata pointer to an MTgomory_ccbk_t structure.

Definition at line 303 of file mt_cplex_cbk.h.

#define MTgomory_ccbk_init __cbdata   ) 
 

Value:

do{\
  MTgomory_ccbk_t*const __MTcbdata = (__cbdata);\
  memset((__MTcbdata),0,sizeof(MTgomory_ccbk_t));\
  MTccbk_info_init(&(__MTcbdata->info));\
  MTsol_init(&(__MTcbdata->sol));\
  MTcut_init(&(__MTcbdata->cut));\
  MTrowmatrix_init(&(__MTcbdata->tb));\
  MTcutHeapInit(&(__MTcbdata->cuth),((size_t)(__MTcbdata->info.max_cpr)));\
  MTlp_init(&(__MTcbdata->lp));}while(0)
initialize an MTgomory_ccbk_t structure.

Parameters:
__cbdata pointer to an MTgomory_ccbk_t structure.

Definition at line 259 of file mt_cplex_cbk.h.

#define MTgomory_ccbk_rsz __cbdata,
__nrowsz,
__ncolsz,
__ntbsz   ) 
 

Value:

do{\
  MTgomory_ccbk_t*const __ccbk = (__cbdata);\
  const int __ccbk_rsz = (__nrowsz);\
  const int __ccbk_csz = (__ncolsz);\
  const int __ccbk_tsz = (__ntbsz);\
  MTrowmatrix_rsz(&(__ccbk->tb),__ccbk_rsz,__ccbk_tsz);\
  MTsol_rsz(&(__ccbk->sol),__ccbk_csz,__ccbk_rsz);\
  MTcut_rsz(&(__ccbk->cut),__ccbk_csz);\
  if(__ccbk->colsz < __ccbk_csz){\
    __ccbk->colsz = __ccbk_csz;\
    EGfree(__ccbk->work);\
    EGfree(__ccbk->extcut);\
    __ccbk->work = EGsMalloc(double,(__ccbk_csz+1));\
    __ccbk->extcut = EGsMalloc(double,__ccbk_csz);}\
  if(__ccbk->rowsz < __ccbk_rsz){\
    __ccbk->rowsz = __ccbk_rsz;\
    EGfree(__ccbk->base);\
    EGfree(__ccbk->f);\
    __ccbk->base = EGsMalloc(double,__ccbk_rsz);\
    __ccbk->f = EGsMalloc(double,__ccbk_rsz);}\
  }while(0)
if required (i.e. __nrowsz > nrowsz || __ncolsz > colsz || __ntbsz > tbsz) , resize an MTgomory_ccbk_t structure

Parameters:
__cbdata pointer to an MTgomory_ccbk_t structure.
__nrowsz new row-size for the structure.
__ncolsz new col-size for the structure.
__ntbsz new tableau-size for the structure.
Note:
When we resize the structure, all previous contents are lost.

Definition at line 277 of file mt_cplex_cbk.h.

#define MTlp_init __lp   )     memset((__lp),0,sizeof(MTlp_t))
 

initialize an MTlp_t structure

Parameters:
__lp pointer to an MTlp_t structure.

Definition at line 106 of file mt_cplex_cbk.h.

#define MTns_ccbk_clear __ns_ccbk   ) 
 

Value:

do{\
  MTns_ccbk_t*const __MTns_ccbk = (__ns_ccbk);\
  EGfree(__MTns_ccbk->rfix);\
  EGfree(__MTns_ccbk->pi);\
  EGfree(__MTns_ccbk->cfix);\
  EGfree(__MTns_ccbk->dj);\
  EGfree(__MTns_ccbk->rstat);\
  EGfree(__MTns_ccbk->sense);\
  EGfree(__MTns_ccbk->oobj);\
  EGfree(__MTns_ccbk->lb);\
  EGfree(__MTns_ccbk->ub);\
  EGfree(__MTns_ccbk->x1);\
  EGfree(__MTns_ccbk->x0);\
  EGfree(__MTns_ccbk->matind);\
  EGfree(__MTns_ccbk->matval);\
  EGfree(__MTns_ccbk->cstat);}while(0)
clear any internally allocated memory within an MTns_ccbk_t structure.

Definition at line 578 of file mt_cplex_cbk.h.

#define MTns_ccbk_display __ns_ccbk   ) 
 

Value:

do{\
  const MTns_ccbk_t*const __MTns = (__ns_ccbk);\
  const double __MTd = ((double)(__MTns->nzrd))/((double)__MTns->ncalls);\
  const double __MTd2 = ((double)(__MTns->loops))/((double)__MTns->ncalls);\
  fprintf(stderr,"\tAverage solution difference  %.7le in %zd calls\n",__MTns->xdiff/((double)__MTns->ncalls),__MTns->ncalls);\
  fprintf(stderr,"\tAverage Primal degenerancy   %.4lf in %zd calls\n",__MTd,__MTns->ncalls);\
  fprintf(stderr,"\tAverage Reoptimization steps %.4lf in %zd calls\n",__MTd2,__MTns->ncalls);}while(0)
Display stats in the given MTns_ccbk_t structure.

Definition at line 596 of file mt_cplex_cbk.h.

#define MTns_ccbk_init __ns_ccbk   ) 
 

Value:

do{\
  MTns_ccbk_t*const __MTns_ccbk = (__ns_ccbk);\
  memset(__MTns_ccbk,0,sizeof(MTns_ccbk_t));}while(0)
initialize an MTns_ccbk_t structure

Definition at line 546 of file mt_cplex_cbk.h.

#define MTns_ccbk_rsz __ns_ccbk,
__ncsz,
__nrsz   ) 
 

Value:

do{\
  MTns_ccbk_t*const __MTns_ccbk = (__ns_ccbk);\
  const size_t __MTns_ncsz = (size_t)(__ncsz);\
  const size_t __MTns_nrsz = (size_t)(__nrsz);\
  if( __MTns_ccbk->rsz < __MTns_nrsz){\
    __MTns_ccbk->rsz = __MTns_nrsz;\
    __MTns_ccbk->rfix = EGrealloc(__MTns_ccbk->rfix,sizeof(char)*__MTns_nrsz);\
    __MTns_ccbk->rstat = EGrealloc(__MTns_ccbk->rstat,sizeof(int)*__MTns_nrsz);\
    __MTns_ccbk->sense = EGrealloc(__MTns_ccbk->sense,sizeof(char)*__MTns_nrsz);\
    __MTns_ccbk->pi = EGrealloc(__MTns_ccbk->pi,sizeof(double)*__MTns_nrsz);}\
  if( __MTns_ccbk->csz < __MTns_ncsz){\
    __MTns_ccbk->csz = __MTns_ncsz;\
    __MTns_ccbk->oobj = EGrealloc(__MTns_ccbk->oobj,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->lb = EGrealloc(__MTns_ccbk->lb,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->ub = EGrealloc(__MTns_ccbk->ub,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->x0 = EGrealloc(__MTns_ccbk->x0,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->x1 = EGrealloc(__MTns_ccbk->x1,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->cfix = EGrealloc(__MTns_ccbk->cfix,sizeof(char)*__MTns_ncsz);\
    __MTns_ccbk->cstat = EGrealloc(__MTns_ccbk->cstat,sizeof(int)*__MTns_ncsz);\
    __MTns_ccbk->matind = EGrealloc(__MTns_ccbk->matind,sizeof(int)*__MTns_ncsz);\
    __MTns_ccbk->matval = EGrealloc(__MTns_ccbk->matval,sizeof(double)*__MTns_ncsz);\
    __MTns_ccbk->dj = EGrealloc(__MTns_ccbk->dj,sizeof(double)*__MTns_ncsz);\
    }\
  }while(0)
resize node-solve arrays

Definition at line 551 of file mt_cplex_cbk.h.

#define MTrowmatrix_init __mt   )     memset((__mt),0,sizeof(MTrowmatrix_t))
 

Initialize an MTrowmatrix_t structure too an empty matrix.

Parameters:
__mt matrix to initialize

Definition at line 133 of file mt_cplex_cbk.h.

#define MTsol_init __sol   )     memset(__sol,0,sizeof(MTsol_t))
 

initialize a MTsol_t structure

Definition at line 150 of file mt_cplex_cbk.h.

#define UPDATE_F_AND_VAL __f,
__cstat,
__bound,
__val   ) 
 

Value:

do{\
  if((__cstat) == CPX_AT_LOWER){\
    (__f) -= (__val) * (__bound);}\
  else if((__cstat) == CPX_AT_UPPER){\
    (__f) -= (__val) * (__bound);\
    (__val) = -1*(__val);}\
  else{TESTG((rval=1),CLEANUP,"Imposible!, cstat is %d", (__cstat));}}while(0)
simple handler to update f according to the complementation of variables

Definition at line 182 of file mt_cplex_cbk.c.


Enumeration Type Documentation

enum MTcutSelection_t
 

types of cut score functions

Enumerator:
MTcsRAW  use raw cut score
MTcsVIOL  cut violation
MTcsINVLINF  $\displaystyle rhs/\|a\|_{\inf} $ for $ax\geq rhs$ inequalities, and as $\displaystyle -rhs/\|a\|_{\inf}$ for $ax\leq rhs$ inequalities.
MTcsINVL2  $\displaystyle rhs^2/\|a\|_2 $ for $ax\geq rhs$ inequalities, and as $\displaystyle -rhs^2/\|a\|_2$ for $ax\leq rhs$ inequalities.
MTcsDefault 

Definition at line 159 of file mt_types.h.


Function Documentation

int main int  argc,
char **  argv
 

this function reads a problem file, set the cut-callback, and call the MIP solver

Definition at line 250 of file mt_cplex.ex.c.

Here is the call graph for this function:

static int mem_limits void   )  [static]
 

function to handle resource usage limits

Definition at line 80 of file mt_cplex.ex.c.

Here is the call graph for this function:

void MTccbk_info_display const MTccbk_info_t *const   info,
FILE *  out
 

display callback info

Definition at line 1110 of file mt_cplex_cbk.c.

int MTccbk_info_process MTccbk_info_t *const   info,
CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
int *const   action
 

test if we should proceed with the cuting process, and update information

Parameters:
info calllback info to use.
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
action pointer to an integer where we return action, zero to continue cutting, non-zero to stop.
Returns:
zero on succsess, non-zero otherwise.

Definition at line 1116 of file mt_cplex_cbk.c.

int MTcheckTabrow const MTlp_t *const   lp,
const double *  tableau,
int *const   status,
double *const   ratio,
int *const   etnz
 

given a tableau in extended notation, verify that it is indeed a numerically stable and correct tableau.

Parameters:
tableau table to check
lp Internal LP description
status on return zero fi everything is OK, non-zero otherwise.
ratio if status==0, store the ratio for the tableau.
etnz if status == 0, store the number of non-zeros in the tableau
Returns:
0 on success, non-zero otherwise

Definition at line 34 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTcompressTabRow const MTlp_t *const   lp,
const MTsol_t *const   sol,
double *  rowval,
int *  rowind,
double *const   rhs,
int *const   nz
 

given a tableau row in extended format, store it as a compresed form, it include setting row to row - integer_part[row] and updating the rhs.

Parameters:
lp internal LP description,
sol current optimal solution.
rowval on input extended row description, on output, non-zeros coefficients
rowind array of length ncols where we will store the indices of non-zero coefficients in the row
nz on output, non-zeros found
Returns:
zero on success, non-zero otherwise

Definition at line 190 of file mt_cplex_cbk.c.

void MTcompute_f const MTlp_t *const   lp,
const double *const   tableau,
const double *const   x,
const double *const   slack,
double *const   f
 

compute ax = f for a tableau

Parameters:
lp pointer to an MTlp_t structure.
tableau where to write the tableau (length lp->nrows +
x structural variables values
slack logical variables values
f where to store ax = f

Definition at line 402 of file mt_cplex_cbk.c.

int MTcompute_fractional_vars const MTlp_t *const   lp,
const double *const   x,
int *const   fracvars,
int *const   nfrac,
const double  minitgap
 

get the indices of fractional (structural) variables

Parameters:
lp pointer to an MTlp_t structure
x current LP x-value
fracvars array of length (at least) MTlp_t::nrows
nfrac pointer to an integer where the number of fractional variables is returned.
minitgap minimum distance to closest integer for a variable to be considered fractional
Returns:
zero on success, non-zero otherwise

Definition at line 339 of file mt_cplex_cbk.c.

int MTcompute_integer_vars const MTlp_t *const   lp,
const double *const   x,
int *const   intvars,
int *const   nint,
const double  intgap
 

get the indices of (almost) integer variables

Parameters:
lp pointer to an MTlp_t structure
x current LP x-value
intvars array of length (at least) MTlp_t::nrows
nint pointer to an integer where the number of fractional variables is returned.
intgap integrality gap allowd for a variable to be considered integer
Returns:
zero on success, non-zero otherwise

Definition at line 309 of file mt_cplex_cbk.c.

int MTcplex_binv_to_tableau const MTlp_t *const   lp,
const double *const   binvrow,
double *const   tableau
 

given an LP desvription, and CPLEX binvrow, compute the actual tableau (including logicals, which are in position i+ncols of the tableau for the i-th logical/slack variable, but we force that the coefficient of equality logicals is zero)

Parameters:
lp pointer to an MTlp_t structure.
binvrow as returned by cplex (length = lp->nrows);
tableau where to write the tableau (length lp->nrows + lp->ncols).
Returns:
zero on success, non-zero otherwise

Definition at line 424 of file mt_cplex_cbk.c.

void MTcplex_display_compress_tableau FILE *  file,
const MTlp_t *const   lp,
const double *const   matval,
const int *const   matind,
const int  nz
 

display a compress tableau into the given file

Parameters:
file where to write the tableau.
lp all LP-related information (including length of the tableau)
matval values of the compresed tableau
matind indices of the non-zeros in the tableau
length of the compresed tableau
Note:
The function support compressed tableau with logical variables

Also, the function would prefix a ~ for all variables that are at its upper bound.

Definition at line 1230 of file mt_cplex_cbk.c.

void MTcplex_display_tableau FILE *  file,
const MTlp_t *const   lp,
const double *const   tableau
 

display a tableau into the given file

Parameters:
file where to write the tableau.
lp all LP-related information (including length of the tableau)
tableau actual values of the tableau.

Definition at line 1199 of file mt_cplex_cbk.c.

int MTcplex_test_cut const MTcut_t *const   cut,
const MTlp_t *const   lp
 

test if a given cut is valid to the given cplex LP

Parameters:
cut cut to be tested
lp pointer to the asociated MTlp_t structure.
Returns:
zero on success, non-zero otherwise.

Definition at line 1748 of file mt_cplex_cbk.c.

int MTcut_check_dominance const MTcut_t *const   cut1,
const MTcut_t *const   cut2
 

check for dominance between cuts.

Parameters:
cut1 First cut to compare
cut2 Second cut to compare
Returns:
zero if neither dominates the other. 1 if cut1 dominates cut2, -1 if cut2 dominates cut1.

Definition at line 193 of file mt_cutpool.c.

Here is the call graph for this function:

void MTcut_clear MTcut_t *const   MTcut  ) 
 

free any internal memory asociated with an MTcut_t structure

Parameters:
MTcut pointer to an MTcut_t structure

Definition at line 392 of file mt_cutpool.c.

void MTcut_copy MTcut_t *const   dest,
const MTcut_t *const   src
 

copy a cut from one structure to the next

Parameters:
orig cut to copy
dest where to copy the cut, all previous information will be lost.

Definition at line 269 of file mt_cutpool.c.

Here is the call graph for this function:

void MTcut_rsz MTcut_t *const   MTcut,
const int  MTnsz
 

resize a cut structure

Parameters:
MTcut pointer to an MTcut_t structure
MTnsz new space for the cut

Definition at line 380 of file mt_cutpool.c.

int MTcutHeap_add_cut const MTlp_t *const   lp,
MTcutHeap_t *const   cuth
 

add the cut (obtained by MTcutHeap_get_free_cut) to the structure, and possibly discard a worst cut from the heap

Parameters:
cuth structure where we perform the operation
Returns:
zero on success, non-zero otherwise.

Definition at line 288 of file mt_cutpool.c.

Here is the call graph for this function:

MTcut_t * MTcutHeap_pop_cut MTcutHeap_t *const   cuth  ) 
 

pop the worst remaining cut within an MTcutHeap_t structure. Note that the memory asociated with the returned MTcut_t structure is handled by the MTcutHeap_t structure, and should not be freed by the user.

Parameters:
cuth structure where we perform the operation
Returns:
pointer to the worst MTcut_t structure (acording to the socre function selected) still in the heap, NULL if no more cuts are left.

Definition at line 364 of file mt_cutpool.c.

void MTcutHeapClear MTcutHeap_t *const   cuth  ) 
 

Clear any internal memory asociated with an MTcutHeap_t structure. Note that after a call to this function a new call to MTcutHeapInit should be made before using any other function over the strucutre.

Parameters:
cuth structure to clear

Definition at line 52 of file mt_cutpool.c.

static int MTcutHeapCn_compute_val MTcutHeap_t *const   cuth,
MTcutHeapCn_t *const   chcn
[static]
 

given a MTcutHeapCn_t compute the value asociated with the stored cut and store it in MTcutHeapCn_t::hcn, inside dbl_EGeHeapCn_t::val, the function used to evaluate depend on the configuration for the cut heap. The smallest the value, the worst the cut is.

Parameters:
chcn structure where we operate
cuth cutheap controling the cut heap connector
Returns:
zero on success, non-zero otherwise

Definition at line 71 of file mt_cutpool.c.

void MTcutHeapInit MTcutHeap_t *const   cuth,
size_t  sz
 

initialize an MTcutHeap_t structure

Parameters:
cuth structure to initialize
sz how many cuts to store

Definition at line 30 of file mt_cutpool.c.

void MTcutNzAlloc void   ) 
 

Definition at line 162 of file mt_cutpool.c.

void MTcutNzFree void   ) 
 

Definition at line 168 of file mt_cutpool.c.

void MTdisplay_lp const MTlp_t *const   lp,
const MTsol_t *const   sol,
FILE *  stream
 

display an LP to the given file

Parameters:
lp pointer to an MTlp_t structure
sol asociated solution to the LP
stream file where to write the LP

Definition at line 369 of file mt_cplex_cbk.c.

int MTget_best_k_integer_tableau_rows const MTlp_t *const   lp,
const MTsol_t *const   sol,
const int  max_tableau,
CPXCENVptr  env,
CPXLPptr  CPXlp,
MTrowmatrix_t *const   tb,
double *const   f
 

Get the best K tableau (asociated to integer basic variables), note that we check for wrongly computed tableau. Note that we compute Z = f + 0.5+Ax, where A is the tableau and f the fractional value (displaced by 1/2), moreover, all basic variables have zero coefficient in the resulting tableau, and all variables are complemented to their bounds.

Parameters:
lp pointer to an MTlp_t structure
sol asociated node solution tb pointer to an MTrowmap_t structure
f fractional value for each tableau.
max_tableau maximum number of tableaus to consider. used to select the set of best integer tableau rows.
env Cplex environment pointer.
CPXlp cplex lp pointer.
Returns:
zero on success, non-zero otherwise.

Definition at line 1265 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTget_best_k_tableau_rows const MTlp_t *const   lp,
const MTsol_t *const   sol,
const int  max_tableau,
CPXCENVptr  env,
CPXLPptr  CPXlp,
MTrowmatrix_t *const   tb,
double *const   f
 

Get the best K tableau (asociated to fractional variables), note that we check for wrongly computed tableau. Note that we compute Z = f + 0.5+Ax, where A is the tableau and f the fractional value (displaced by 1/2), moreover, all basic variables have zero coefficient in the resulting tableau, and all variables are complemented to their bounds.

Parameters:
lp pointer to an MTlp_t structure
sol asociated node solution tb pointer to an MTrowmap_t structure
f fractional value for each tableau.
max_tableau maximum number of tableaus to consider. used to select the set of best fractional tableau rows.
env Cplex environment pointer.
CPXlp cplex lp pointer.
Returns:
zero on success, non-zero otherwise.

Definition at line 1354 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTget_sol const MTlp_t *const   lp,
CPXCENVptr  env,
CPXLPptr  CPXlp,
void *  cbdata,
int  wherefrom,
MTsol_t *const   sol
 

obtain a solution to the current node LP

Parameters:
lp pointer to an MTlp_t structure
env Cplex environment pointer.
CPXlp cplex lp pointer.
sol where to store the solution
CPLEX callback cbdata.
wherefrom CPLEX callback wherefrom parameter.
Returns:
zero on success, non-zero otherwise

Definition at line 1714 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTgomory_ccbk CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
void *  cbhandle,
int *  useraction_p
 

Interface for the CPLEX-callback for the MTgomoryCut function.

Parameters:
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
cbhandle A pointer to user private data.
useraction_p A pointer to an integer indicating the action for ILOG CPLEX to take at the completion of the user callback. The table summarizes possible actions.
  • 0 CPX_CALLBACK_DEFAULT Use cuts as added
  • 1 CPX_CALLBACK_FAIL Exit optimization
  • 2 CPX_CALLBACK_SET Use cuts as added

Definition at line 579 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTgomoryCut const int  nrows,
const int  ncols,
const double *const   rowval,
const int *const   rowind,
const int *const   rowbeg,
const double *const   f,
const double *const   base,
double *const   cutval,
double *const   work
 

Given a set of tableau rows in column format, return the asociated multi-tableau gomory cut. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols. The returned constraint is $\displaystyle \sum\limits_{i\in I}s_ib_i\geq 1$ where $ \displaystyle b_i=$ cutval[i], $\displaystyle \forall i\in I$. And where the set used to cut is $\displaystyle L=\left\{x:c\sigma x\leq c_o ,\forall\sigma\in\mathrm{diag}\left(\{-1,1\}^{m}\right)\right\}$ and $\displaystyle c=$ base, $\displaystyle c_o = \frac12\sum\limits_{i=1}^mc_i$, where we assume that $c_i\neq 0\forall i=1,\ldots,m$. Finally, note that $\displaystyle b_i=\max\limits_{\sigma:c\sigma r_i>0}\left\{\frac{c\sigma r_i}{c_o-c\sigma f}\right\}$.

Parameters:
nrows number of rows in the tableau description
ncols number of variables in the tableau description
rowval array of values of the tableau.
rowind array of column-indices for the corresponding value in rowval.
rowbeg array of length nrows+1 storing the position in which position of rowval the j-th row non-zero coefficient beggins, and in position nrows it is stored the number of non-zeros in the description.
f fractional value asociated to the tableau.
cutval array of length ncols storing the coefficient for each variable in the cut (the array is assumed to have length at least ncols).
base array of length nrows, where each part coefficient is assumed to be non-zero.
work array of length at least ncols+1 to be used as working space.
Returns:
zero on succes, non-zero otherwise.

Definition at line 33 of file mt_gomory.c.

Here is the call graph for this function:

void MTintPermSort const int  sz,
int *const   perm,
const int *const   elem
 

permutation sort for integers

Definition at line 129 of file mt_cutpool.c.

void MTlp_clear MTlp_t *const   MTlp  ) 
 

free any internally allocated memory inside an MTlp_t structure.

Parameters:
MTlp pointer to an MTlp_t structure.

Definition at line 1175 of file mt_cplex_cbk.c.

int MTlp_load_problem CPXCENVptr  env,
MTlp_t *const   lp,
CPXLPptr  CPXlp,
void *  cbdata,
int  wherefrom
 

Given an LP problem pointer, load it into an (initialized) MTlp_t structure.

Parameters:
lp pointer to an MTlp_t structure.
CPXlp pointer to a CPLEX LP problem
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
Returns:
zero on success.

Definition at line 469 of file mt_cplex_cbk.c.

Here is the call graph for this function:

static void MTlp_rsz MTlp_t *const   MTlp,
const int  rsz,
const int  csz,
const int  tsz
[inline, static]
 

resize an lp structure

Parameters:
MTlp pointer to an MTlp_t structure.
MTrsz new row-space (it will only grow).
MTcsz new col space (it will only grow).
MTtsz new non-zero space (it will only grow).

Definition at line 259 of file mt_cplex_cbk.c.

int MTnode_solve CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
void *  cbhandle,
int *  useraction_p
 

simple node solve callback, we use it to count primal degenerancy, i.e. the number of variables with zero reduced cost, this should indicate posibility of doing nested objective funcions to avoid primal degenerancy

Definition at line 1848 of file mt_cplex_cbk.c.

int MTraw_to_complemented_cut const MTlp_t *const   lp,
const MTsol_t *const   sol,
double *const   extended,
MTcut_t *const   cut
 

transform an extened non-complemented cut to a compresed and complemented form. We assume the extended cut is of the form ax >= 1. Note that we expect no coefficient for basic variables in the cut, moreover, the rhs of the inequality is assumed to be 1 (in raw format).

Note:
The parameter extended WILL CHANGE ITS CONTENTS
Parameters:
cut pointer to an MTcut_t structure, we will store the compresed cut here.
sol asociated node solution
extended extended cut (length MTlp_t::nrows + MTlp_t::ncols), note that this informnation will be modified, and will contain the complemented cut in extended form.
lp pointer to the related MTlp_t structure
Returns:
zero on success, non-zero otherwise.

Definition at line 1553 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTraw_to_uncomplemented_cut const MTlp_t *const   lp,
const double *const   extended,
const int  nactive,
const int *const   active,
MTcut_t *const   cut,
int *const   status
 

transform an extened non-complemented cut to a compresed for. We assume the extended cut is of the form ax >= 1.m Note that we expect no coefficient for basic variables in the cut.

Note:
The parameter extended WILL NOT CHANGE ITS CONTENTS
Parameters:
cut pointer to an MTcut_t structure, we will store the compresed cut here.
extended extended cut (length MTlp_t::nrows + MTlp_t::ncols), note that this informnation will be modified, and will contain the complemented cut in extended form.
nactive length of array active
active if not null, it should store the positions of possibly active positions in extended; all positions outside it are assumed to be zero.
lp pointer to the related MTlp_t structure
status if status is zero, the cut should be stored, otherwise, it should be discarded.
Returns:
zero on success, non-zero otherwise.

Definition at line 1443 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTread_cplex_options CPXENVptr  env,
FILE *  inputfile
 

read optionfile and set parameter to CPLEX

Parameters:
env pointer to the CPLEX environment structure.
inputfile pointer to the optionfile.

Definition at line 701 of file mt_cplex_cbk.c.

int MTreadRTableau FILE *  input,
int *const   nrows,
int *const   ncols,
double **const   rowval,
int **const   rowind,
int **const   rowbeg,
double **const   f
 

read a tableau in row format from a file (see MTgomoryCut for details on the tableau information).

Parameters:
input file where to write the tableau
nrows number of rows in the tableau description
ncols number of variables in the tableau description
rowval array of values of the tableau.
rowind array of column-indices for the corresponding value in rowval.
rowbeg array of length nrows+1 storing the position in which position of rowval the j-th row non-zero coefficient beggins, and in position nrows it is stored the number of non-zeros in the description.
f fractional value asociated to the tableau.
Note:
We assume that the pointers given contain no meaningfull data

Definition at line 59 of file mt_tableau.c.

void MTrowmatrix_clear MTrowmatrix_t *const   MTrm  ) 
 

free any internally allocated memory inside an MTrowmatrix_t structure.

Parameters:
MTrm rowmatrix to clear

Definition at line 1166 of file mt_cplex_cbk.c.

void MTrowmatrix_rsz MTrowmatrix_t *const   MTrm,
const int  MTrsz,
const int  MTnsz
 

resize a MTrowmatrix_t structure

Parameters:
MTrm rowmatrix to resize
MTrsz new row-space (it will only grow).
MTnsz new non-zero space (it will only grow).

Definition at line 1092 of file mt_cplex_cbk.c.

int MTsl1 const double *const   ineq,
MTrowmatrix_t *const   tb,
double *const   work,
const char *const   vtype,
const int  mode
 

compute in work the value ineq * r_i as described by the given inequality and tableau, we assume that work is zero, moreover, given the variable status (interger/binary) and an operation mode, try to do nothing (mode 0), try to get ineq*r[i] >0 (mode 1) or try to get ineq*r[i] < 0 (mode 2), note that mode != 0 will modify the tableau

Note:
We assume that the coefficients for integer variables in the tableau are between ]-1:1[

when computing a new inequality, mode can be non-zero only at the beggining, otherwise, the resulting inequality might be invalid.

Returns:
zero on success, non-zero otherwise

Definition at line 41 of file mt_SL.c.

int MTsl2 const double *const   ineq,
double *const   cutval,
double *const   work,
const int *const   active,
const int  nactive
 

, set cutval[i] = max(cutval[i],work[i] > 0 ? work[i]/ineq[3]: 0); and reset work to zero

Definition at line 90 of file mt_SL.c.

int MTsl3 const double *const   f,
double *const   ineq,
const int  dim
 

given an inequality and f compute ineq[dim+1] and (possible change signs) to ensure that ineq[dim+1] > 0;

Definition at line 120 of file mt_SL.c.

int MTsl4 double *const   ineq,
const double  x0,
const double  _y0,
const double  x1,
const double  _y1
 

compute the line that pases through the two given points

Definition at line 142 of file mt_SL.c.

int MTsl5 int *const   nactive,
int *const   active,
const int  nvars,
const int  nz,
const double *const   compval,
const int *const   compind
 

compute active set from the given compresed vector set

Definition at line 160 of file mt_SL.c.

int MTsl_ccbk CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
void *  cbhandle,
int *  useraction_p
 

Interface for the CPLEX-callback for the MTt1Cut function.

Parameters:
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
cbhandle A pointer to user private data.
useraction_p A pointer to an integer indicating the action for ILOG CPLEX to take at the completion of the user callback. The table summarizes possible actions.
  • 0 CPX_CALLBACK_DEFAULT Use cuts as added
  • 1 CPX_CALLBACK_FAIL Exit optimization
  • 2 CPX_CALLBACK_SET Use cuts as added

Definition at line 458 of file mt_SL.c.

Here is the call graph for this function:

int MTslCut MTgomory_ccbk_t *const   data,
CPXCENVptr  env,
CPXLPptr  lp
 

The idea of these cuts is to use some ideas from the paper ``On the relative strengt of split, triangle and cuadrilateral cuts'' from Cornuejols et. al. The main point is that splits are bad (compared against triangles or quadrilaterals) whenever we have (an almost) integer part in $f_i$ for some coordiante $i$, and a wide base, this can be seen on the following figure: $ \hspace*{\textwidth} \vspace*{2in} x_1 = 1 \psset{arrows=->,fillstyle=solid,fillcolor=black,linecolor=black} \begin{pspicture}(0,0)(5,3) \multido{\i=0+1}{6}{{ \pscircle(\i,0){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,1){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,2){1mm} }} \multido{\i=0+1}{6}{{ \pscircle(\i,3){1mm} }} \rput(2.5,2){$f$} \psline(2.5,2)(0,.5)\nbput{$r_1$} \psline(2.5,2)(4,1)\nbput{$r_2$} \psline[arrows=|-|](1,1)(4,1)\nbput{M} \end{pspicture} $.

Parameters:
data calback data, we assume that the solution and lp has already been loaded into the structure.
env CPLEX environment
lp CPLEX lp pointer to the current LP
Returns:
zero on success, non-zero otherwise.

Definition at line 177 of file mt_SL.c.

Here is the call graph for this function:

void MTsol_clear MTsol_t *const   MTsol  ) 
 

clear any internal memory asociated with an MTsol_t structure

Parameters:
MTsol structure to clear

Definition at line 1157 of file mt_cplex_cbk.c.

void MTsol_rsz MTsol_t *const   MTsol,
const int  cnz,
const int  rnz
 

resize a solution structure

Parameters:
MTsol solution to resize
cnz new column space
rnz new row space

Definition at line 1074 of file mt_cplex_cbk.c.

int MTt1_ccbk CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
void *  cbhandle,
int *  useraction_p
 

Interface for the CPLEX-callback for the MTt1Cut function.

Parameters:
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
cbhandle A pointer to user private data.
useraction_p A pointer to an integer indicating the action for ILOG CPLEX to take at the completion of the user callback. The table summarizes possible actions.
  • 0 CPX_CALLBACK_DEFAULT Use cuts as added
  • 1 CPX_CALLBACK_FAIL Exit optimization
  • 2 CPX_CALLBACK_SET Use cuts as added

Definition at line 169 of file mt_T1.c.

Here is the call graph for this function:

int MTt1Cut const MTrowmatrix_t *const   tb,
MTlp_t *const   lp,
MTcutHeap_t *const   cuth,
const double *const   f,
double *const   extcut,
double *const   work,
const int  cut_style
 

compute the best cuts found using T1 sets as a base, where

\[ T1_n=\left\{ x\in\mathbb{R}^n: x_i\geq 0, \sum\limits_{i=1}^n x_i\leq n \right\},\]

and where we try all possible $2^n$ reorientation of the axis. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols.

Parameters:
tb tableau description to use.
lp description of the current node LP.
f fractional values asociated with the tableau.
work array of length at least ncols+1 to be used as working space.
cuth cutheap to store the best cuts found in the process
extcut array of length ncols to store extended cut while processing
cut_style if 0, produce cut for cuts derived with exactly tb->nrows rows, if 1, also try with smaller values up to 2.
Returns:
zero on success, non-zero otherwise.

Definition at line 30 of file mt_T1.c.

Here is the call graph for this function:

int MTt2_ccbk CPXCENVptr  env,
void *  cbdata,
int  wherefrom,
void *  cbhandle,
int *  useraction_p
 

Interface for the CPLEX-callback for the MTt2Cut function.

Parameters:
env as returned by CPXopenCPLEX function.
cbdata A pointer passed from the optimization routine to the user-written callback that identifies the problem being optimized. The only purpose of this pointer is to pass it to the callback information routines.
wherefrom An integer value indicating where in the optimization this function was called. It has the value CPX_CALLBACK_MIP_CUT.
cbhandle A pointer to user private data.
useraction_p A pointer to an integer indicating the action for ILOG CPLEX to take at the completion of the user callback. The table summarizes possible actions.
  • 0 CPX_CALLBACK_DEFAULT Use cuts as added
  • 1 CPX_CALLBACK_FAIL Exit optimization
  • 2 CPX_CALLBACK_SET Use cuts as added

Definition at line 199 of file mt_T2.c.

Here is the call graph for this function:

int MTt2Cut const MTrowmatrix_t *const   tb,
MTlp_t *const   lp,
MTcutHeap_t *const   cuth,
const double *const   f,
double *const   extcut,
double *const   work,
const int  cut_style
 

compute the best cuts found using T2 sets as a base, where

\[ T2_n=\left\{x\in\mathbb{R}^n: \begin{matrix} R1:&-x_1 & \leq 0\\ R2:& x_1-x_2 & \leq 1\\ R3:& x_1+x_2-x_3 & \leq 2\\ \vdots & \vdots & \vdots\\ Rn:& x_1+\ldots+x_{n-1}-x_n&\leq n-1\\ R(n+1):& x_1+\ldots+x_{n-1}+x_n&\leq n\\ \end{matrix}\right\} \]

, and where we try all possible $n!2^n$ reorientation of the axis and permutation of the rows. Note that we assume that all variables are continuous, and that the tableau represent the (vectorial) equation

\[\vec{x} = \vec{f}+\frac12\vec{1} +\sum\limits_{i\in I}s_i\vec{r_i}, \]

where $\displaystyle\vec{x}\in\mathbb{Z}^m$, $\displaystyle\vec{f}\in\left]-\frac12,\frac12\right[^m$, $\displaystyle\vec{r_i}\in\mathbb{R}^m,\forall i\in I$, $s\displaystyle _i\in\mathbb{R}_+,\forall i\in I$; and where $\displaystyle m=$ nrows and $\displaystyle |I|=$ ncols.

Parameters:
tb tableau description to use.
lp description of the current node LP.
f fractional values asociated with the tableau.
work array of length at least ncols+1 to be used as working space.
cuth cutheap to store the best cuts found in the process
extcut array of length ncols to store extended cut while processing
cut_style if 0, produce cut for cuts derived with exactly tb->nrows rows, if 1, also try with smaller values up to 2.
Returns:
zero on success, non-zero otherwise.

Definition at line 30 of file mt_T2.c.

Here is the call graph for this function:

int MTuncomplemented_to_complemented_cut const MTlp_t *const   lp,
const MTsol_t *const   sol,
MTcut_t *const   cut
 

transform an uncomplemented cut to a complemented form. Note that we expect no coefficient for basic variables in the cut.

Parameters:
cut pointer to an MTcut_t structure, we will complement the cut.
sol asociated node solution
lp pointer to the related MTlp_t structure
Returns:
zero on success, non-zero otherwise.

Definition at line 1532 of file mt_cplex_cbk.c.

Here is the call graph for this function:

int MTwriteRTableau FILE *  out,
const int  nrows,
const int  ncols,
const double *const   rowval,
const int *const   rowind,
const int *const   rowbeg,
const double *const   f
 

write a tableau in row format into a file (see MTgomoryCut for details on the tableau information).

Parameters:
out file where to write the tableau
nrows number of rows in the tableau description
ncols number of variables in the tableau description
rowval array of values of the tableau.
rowind array of column-indices for the corresponding value in rowval.
rowbeg array of length nrows+1 storing the position in which position of rowval the j-th row non-zero coefficient beggins, and in position nrows it is stored the number of non-zeros in the description.
f fractional value asociated to the tableau.

Definition at line 32 of file mt_tableau.c.

int parseargs int  argc,
char **  argv
 

parsing parameters function, also read problem and options file

Definition at line 163 of file mt_cplex.ex.c.

Here is the call graph for this function:

static void sighandler int  s  )  [static]
 

Definition at line 66 of file mt_cplex.ex.c.

void usage char *const   prog  ) 
 

usage function

Definition at line 127 of file mt_cplex.ex.c.


Variable Documentation

char __mt_cplex_errbuf[4096]
 

internal error string for cplex

Definition at line 400 of file mt_cplex_cbk.c.

char __mt_cplex_errbuf[4096]
 

internal error string for cplex

Definition at line 400 of file mt_cplex_cbk.c.

int cut_dominance = 0 [static]
 

cut dominance check enable

Definition at line 51 of file mt_cplex.ex.c.

MTcutSelection_t cut_selection = MTcsDefault [static]
 

cut selection strategy

Definition at line 49 of file mt_cplex.ex.c.

int cut_style = 1 [static]
 

cut generation strategy

Definition at line 59 of file mt_cplex.ex.c.

int cuts_at_root_only = 1 [static]
 

cuts at root only

Definition at line 47 of file mt_cplex.ex.c.

CPXENVptr env = 0 [static]
 

global cplex environment

Definition at line 39 of file mt_cplex.ex.c.

CPXLPptr lp = 0 [static]
 

global lp pointer

Definition at line 41 of file mt_cplex.ex.c.

int max_cuts_per_node = 1000 [static]
 

max cuts per node

Definition at line 43 of file mt_cplex.ex.c.

int max_cuts_per_round = 10 [static]
 

max cuts per round

Definition at line 45 of file mt_cplex.ex.c.

double max_rtime = INT_MAX [static]
 

maximum running time

Definition at line 61 of file mt_cplex.ex.c.

int max_tableau_rows = 5 [static]
 

maximum number of tableau-rows to use

Definition at line 53 of file mt_cplex.ex.c.

unsigned long memlimit = UINT_MAX [static]
 

maximum memory usage

Definition at line 63 of file mt_cplex.ex.c.

double MTccbk_discarded_ratio
 

minimum ratio of discarded inequality

Definition at line 577 of file mt_cplex_cbk.c.

double MTccbk_discarded_ratio = DBL_MAX
 

minimum ratio of discarded inequality

Definition at line 577 of file mt_cplex_cbk.c.

double MTccbk_discarded_violation
 

maximum violation of discarded inequality

Definition at line 576 of file mt_cplex_cbk.c.

double MTccbk_discarded_violation = 0
 

maximum violation of discarded inequality

Definition at line 576 of file mt_cplex_cbk.c.

int MTccbk_fail_ratio
 

Number of discarded cuts because high ratio.

Definition at line 575 of file mt_cplex_cbk.c.

int MTccbk_fail_ratio = 0
 

Number of discarded cuts because high ratio.

Definition at line 575 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_btype
 

Type of basic variable for tableau is not integer.

Definition at line 572 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_btype = 0
 

Type of basic variable for tableau is not integer.

Definition at line 572 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_coeff
 

Number of discarded cuts because tableaus has wrong coefficient for basic variable.

Definition at line 570 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_coeff = 0
 

Number of discarded cuts because tableaus has wrong coefficient for basic variable.

Definition at line 570 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_frac
 

Number of discarded cuts because fraction is too small.

Definition at line 571 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_frac = 0
 

Number of discarded cuts because fraction is too small.

Definition at line 571 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_nbasic
 

Number of discarded cuts because tableaus has wrong number of basic variables.

Definition at line 573 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_nbasic = 0
 

Number of discarded cuts because tableaus has wrong number of basic variables.

Definition at line 573 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_ratio
 

display summary information for the call-back calls

Parameters:
__fout where to write the sumary

Definition at line 569 of file mt_cplex_cbk.c.

int MTccbk_fail_tableau_ratio = 0
 

Definition at line 569 of file mt_cplex_cbk.c.

int MTccbk_fail_violation
 

Number of discarded cuts because low violation.

Definition at line 574 of file mt_cplex_cbk.c.

int MTccbk_fail_violation = 0
 

Number of discarded cuts because low violation.

Definition at line 574 of file mt_cplex_cbk.c.

int* MTcut1_nz = 0 [static]
 

internal memory for dominance check

Definition at line 160 of file mt_cutpool.c.

int* MTcut2_nz = 0 [static]
 

Definition at line 161 of file mt_cutpool.c.

const char MTversion[1024]
 

Version String of the program.

char options[1024] = "cplex.opt" [static]
 

CPLEX option file name.

Definition at line 35 of file mt_cplex.ex.c.

char problem[1024] = "" [static]
 

input problem file

Definition at line 37 of file mt_cplex.ex.c.

int use_cuts = 1 [static]
 

use callback cuts

Definition at line 55 of file mt_cplex.ex.c.

int use_solve_cbk = 0 [static]
 

use solve callback

Definition at line 57 of file mt_cplex.ex.c.


Generated on Mon Oct 26 09:17:23 2009 for MTgomory by  doxygen 1.4.6