#include <stdlib.h>
#include <string.h>
#include <float.h>
#include <math.h>
#include <assert.h>
#include <stdio.h>
#define SMUT_NO_IO
#include "smut.h"
#undef SMUT_NO_IO
#define AUCTION_INTERNAL
#include "auction.h"
#undef AUCTION_INTERNAL
#if !defined(HAVE_RESTRICT)
#define restrict
#endif
#if !defined(HAVE_INLINE)
#define inline
#endif
#if !defined(NOSTATS)
#define STATS(x) do { if (a_nstats < A_MAX_STATS) x } while (0)
#define A_MAX_STATS 1024
int a_nstats;
double a_stats_mu_storage[A_MAX_STATS];
double *a_stats_mu = &a_stats_mu_storage[0];
int a_stats_n_unmatched_storage[A_MAX_STATS];
int *a_stats_n_unmatched = &a_stats_n_unmatched_storage[0];
double a_stats_gap_storage[A_MAX_STATS];
double *a_stats_gap = &a_stats_gap_storage[0];
#else
#define STATS(x) do { } while (0)
#endif
void
auction_toexpint (const struct spcsr_t *A, double *expent_in)
{
int k;
const int nent = A->nent;
const double * restrict entry = A->entry;
double * restrict expent = expent_in;
for (k = 0; k < nent; ++k) {
const double e = entry[k];
#if defined(HAVE_LOGB)
expent[k] = logb(fabs(e));
#elif defined(HAVE_LOG2)
expent[k] = floor(log2(fabs(e)));
#else
{
int exp;
frexp(fabs(e), &exp);
expent[k] = exp;
}
#endif
}
}
void
auction_toexp (const struct spcsr_t *A, double *expent_in)
{
int k;
const int nent = A->nent;
const double * restrict entry = A->entry;
double * restrict expent = expent_in;
for (k = 0; k < nent; ++k) {
const double e = entry[k];
expent[k] = log2(fabs(e));
}
}
void
auction_tounit (const struct spcsr_t *A, double *unit)
{
int k;
const int nent = A->nent;
const double * restrict entry = A->entry;
for (k = 0; k < nent; ++k) {
if (entry[k] != 0.0)
unit[k] = 1.0;
else
unit[k] = 0.0;
}
}
void
auction_tounitall (const struct spcsr_t *A, double *unit)
{
int k;
const int nent = A->nent;
for (k = 0; k < nent; ++k)
unit[k] = 1.0;
}
void
auction_toabs (const struct spcsr_t *A, double *absent_in)
{
int k;
const int nent = A->nent;
const double * restrict entry = A->entry;
double * restrict absent = absent_in;
for (k = 0; k < nent; ++k)
absent[k] = fabs(entry[k]);
}
void
auction_toabsint (const struct spcsr_t *A, double *absent_in)
{
int k;
const int nent = A->nent;
const double * restrict entry = A->entry;
double * restrict absent = absent_in;
for (k = 0; k < nent; ++k)
absent[k] = floor(fabs(entry[k]));
}
double
auction_shift (const struct spcsr_t *A, double *b_in)
{
int k;
const int nent = A->nent;
double * restrict b = b_in;
/* double mx = -HUGE_VAL; */
double mn = HUGE_VAL;
for (k = 0; k < nent; ++k) {
const double e = b[k];
/*if (e > mx) mx = e;*/
if (e > -HUGE_VAL && e < mn) mn = e;
}
for (k = 0; k < nent; ++k) {
b[k] -= mn;
}
/*fprintf (stderr, "shifting by %d\n", mn);*/
return mn;
}
void
auction_scaleshift (const struct spcsr_t *A, double *b_in)
{
int i, j, k;
const int nr = A->nr;
const int nc = A->nc;
const int nent = A->nent;
const int * restrict rowoff = A->rowoff;
const int * restrict colind = A->colind;
double * restrict b = b_in;
double * restrict Cmin;
Cmin = malloc(nc * sizeof(double));
for (j = 0; j < nc; ++j) Cmin[j] = HUGE_VAL;
for (i = 0; i < nr; ++i) {
double mn = HUGE_VAL;
for (k = rowoff[i]; k < rowoff[i+1]; ++k) {
const double e = b[k];
if (e > -HUGE_VAL && e < mn) mn = e;
}
if (mn == HUGE_VAL) mn = 0.0;
for (k = rowoff[i]; k < rowoff[i+1]; ++k) {
const int j = colind[k];
double e = b[k];
e -= mn;
b[k] = e;
if (e > -HUGE_VAL && e < Cmin[j]) Cmin[j] = e;
}
}
for (j = 0; j < nc; ++j) if (Cmin[j] == HUGE_VAL) Cmin[j] = 0.0;
for (k = 0; k < nent; ++k) {
const int j = colind[k];
b[k] -= Cmin[j];
}
free (Cmin);
}
void
auction_eval_primal_mdual (const struct spcsr_t A, const double *b_in,
const int *match_in, const double *price_in,
double *primal_out, double *mdual_out)
{
const int nr = A.nr;
int * restrict rowoff = A.rowoff;
int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
const int * restrict match = match_in;
double primal = 0.0, mdual = 0.0;
int i;
/*
take care only to add only those prices in the matching.
then is can be used for a parallel routine.
*/
for (i = 0; i < nr; ++i) {
int k;
double profit = -HUGE_VAL;
for (k = rowoff[i]; k < rowoff[i+1]; ++k) {
const int j = colind[k];
const double e = b[k];
const int mi = match[j];
const double p = price[j];
double v;
v = e - p;
if (mi == i) {
primal += e;
if (p < HUGE_VAL)
mdual += p;
else
mdual += e; /* "inf" prices cancel with "-inf" profits */
}
if (v > profit) profit = v;
}
if (profit > -HUGE_VAL) {
mdual += profit;
}
}
*primal_out = primal;
*mdual_out = mdual;
}
void
auction_eval_primal_dual (const struct spcsr_t A, const double *b_in,
const int *match_in, const double *price_in,
double *primal_out, double *dual_out)
{
const int nr = A.nr;
const int nc = A.nc;
int * restrict rowoff = A.rowoff;
int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
const int * restrict match = match_in;
double primal = 0.0, dual = 0.0;
int i, j;
/*
this one adds _all_ the prices into the dual.
*/
for (j = 0; j < nc; ++j) if (price[j] < HUGE_VAL) dual += price[j];
for (i = 0; i < nr; ++i) {
int k;
double profit = -HUGE_VAL;
for (k = rowoff[i]; k < rowoff[i+1]; ++k) {
const int j = colind[k];
const double e = b[k];
const int mi = match[j];
const double p = price[j];
double v;
v = e - p;
if (mi == i) {
primal += e;
if (HUGE_VAL == p)
dual += e; /* "inf" prices cancel with "-inf" profits */
}
if (v > profit) profit = v;
}
if (profit > -HUGE_VAL) {
dual += profit;
}
}
*primal_out = primal;
*dual_out = dual;
}
double
auction_eval_gap (const struct spcsr_t A, const double *b_in,
const int *match_in, const double *price_in,
double *gap_out)
{
const int nr = A.nr;
int * restrict rowoff = A.rowoff;
int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
const int * restrict match = match_in;
double gap = 0.0;
double primal = 0.0;
int i;
for (i = 0; i < nr; ++i) {
int k;
double profit = -HUGE_VAL;
double match_e = HUGE_VAL;
double match_p = HUGE_VAL;
for (k = rowoff[i]; k < rowoff[i+1]; ++k) {
const int j = colind[k];
const double e = b[k];
const int mi = match[j];
const double p = price[j];
double v;
v = e - p;
if (mi == i) {
match_e = e;
match_p = p;
}
if (v > profit) profit = v;
}
if (match_p < HUGE_VAL)
gap += (match_p + profit) - match_e;
if (!isinf(match_e)) primal += fabs(match_e);
}
*gap_out = gap;
return primal;
}
#if 1
inline struct auction_bid
auction_find_bid (const struct spcsr_t A, const double *b_in,
const int i,
const double *price_in)
{
int k;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
struct auction_bid out;
double best_val_0 = -HUGE_VAL;
double best_val_1 = -HUGE_VAL;
double best_ent = -HUGE_VAL;
int best_match = -1;
const int offend = A.rowoff[i+1];
for (k = rowoff[i]; k < offend; ++k) {
const int j = colind[k];
const double e = b[k];
const double p = price[j];
double val;
val = e - p;
#if 0
{
//extern int rank;
if ((2518 == i || 680 == i || j == 2050 || 881 == i)) { // && 0 == rank) {
fprintf (stderr, "%d: %d: j %d, e %g, price %g (%g: %g %g)\n", 0, //rank,
i, j, e, p,
val, best_val_1, best_val_0);
}
}
#endif
if (val <= best_val_1) continue;
if (val <= best_val_0) {
best_val_1 = val;
continue;
}
best_val_1 = best_val_0;
best_val_0 = val;
best_match = j;
best_ent = e;
}
/*
{
extern int rank;
if (best_match == -1) {
fprintf (stderr, "%d: row %d, [%d, %d)\n", rank,
i, rowoff[i], rowoff[i+1]);
}
}
*/
out.price = best_ent - best_val_1;
out.ent = best_ent;
out.val = best_val_0;
out.j = best_match;
return out;
}
#else
inline struct auction_bid
auction_find_bid (const struct spcsr_t A, const double *b_in,
const int i,
const double *price_in)
{
int k;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
struct auction_bid out;
double best_val_0 = -HUGE_VAL;
double best_val_1 = -HUGE_VAL;
double best_ent = -HUGE_VAL;
int best_match = -1;
const int offend = A.rowoff[i+1];
for (k = rowoff[i]; k < offend; ++k) {
const int j = colind[k];
const double e = b[k];
const double p = price[j];
double val;
val = e - p;
if (val > best_val_1) {
if (val >= best_val_0) {
best_val_1 = best_val_0;
best_val_0 = val;
best_match = j;
best_ent = e;
}
else
best_val_1 = val;
}
}
out.price = best_ent - best_val_1;
out.ent = best_ent;
out.val = best_val_0;
out.j = best_match;
return out;
}
#endif
inline struct auction_bid
auction_find_bid_cache (const struct spcsr_t A, const double *b_in,
const int i,
const double *price_in,
int *jcache_in, double *ecache_in)
{
int k;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
const double * restrict b = b_in;
const double * restrict price = price_in;
int * restrict jcache = jcache_in; /* i i / i i / ... */
double * restrict ecache = ecache_in; /* e e v / e e v / ... */
struct auction_bid out;
double best_val_0;
double best_val_1;
double best_val_2;
double best_ent;
double best_ent_1;
int best_match;
int best_match_1;
double p0, p1;
int offend;
best_match = jcache[2*i];
best_match_1 = jcache[2*i+1];
p0 = price[best_match];
p1 = price[best_match_1];
best_ent = ecache[3*i];
best_ent_1 = ecache[3*i+1];
best_val_2 = ecache[3*i+2];
best_val_0 = best_ent - p0;
best_val_1 = best_ent_1 - p1;
if (best_val_0 > best_val_2 && best_val_1 > best_val_2) {
/*
if (581 == i)
fprintf (stderr, "cache hit %d: %g(%g-%g) %g(%g-%g) %g\n", i,
best_val_0, best_ent, p0, best_val_1, best_ent_1, p1, best_val_2);
*/
if (best_val_0 > best_val_1) {
out.price = best_ent - best_val_1;
/*out.price = p0 + best_val_0 - best_val_1;*/
out.ent = best_ent;
out.val = best_val_0;
out.j = best_match;
}
else {
out.price = best_ent_1 - best_val_0;
/*out.price = p1 + best_val_1 - best_val_0;*/
out.ent = best_ent_1;
out.val = best_val_1;
out.j = best_match_1;
}
return out;
}
/*
if (581 == i)
fprintf (stderr, "cache miss %d: %g(%g-%g) %g(%g-%g) %g\n", i,
best_val_0, best_ent, p0, best_val_1, best_ent_1, p1, best_val_2);
*/
offend = A.rowoff[i+1];
best_val_0 = -HUGE_VAL;
best_val_1 = -HUGE_VAL;
best_val_2 = -HUGE_VAL;
best_ent = -HUGE_VAL;
best_ent_1 = -HUGE_VAL;
best_match = -1;
best_match_1 = -1;
for (k = rowoff[i]; k < offend; ++k) {
const int j = colind[k];
const double e = b[k];
const double p = price[j];
double val;
val = e - p;
/*
{
extern int rank;
if (2027 == i && 1 == rank) {
fprintf (stderr, "%d: %d: j %d, e %g, price %g\n", rank,
i, j, e, p);
}
}
*/
if (val <= best_val_2) continue;
if (val <= best_val_1) {
best_val_2 = val;
continue;
}
if (val <= best_val_0) {
best_val_2 = best_val_1;
best_val_1 = val;
best_match_1 = j;
best_ent_1 = e;
continue;
}
best_val_2 = best_val_1;
best_match_1 = best_match;
best_ent_1 = best_ent;
best_val_1 = best_val_0;
best_val_0 = val;
best_match = j;
best_ent = e;
}
/*
{
extern int rank;
if (best_match == -1) {
fprintf (stderr, "%d: row %d, [%d, %d)\n", rank,
i, rowoff[i], rowoff[i+1]);
}
}
*/
/*
if (581 == i)
fprintf (stderr, "cache miss2 %d: %g(%g-%g) %g(%g-%g) %g\n", i,
best_val_0, best_ent, price[best_match],
best_val_1, best_ent_1, price[best_match_1], best_val_2);
*/
jcache[2*i] = best_match;
jcache[2*i+1] = best_match_1;
ecache[3*i] = best_ent;
ecache[3*i+1] = best_ent_1;
ecache[3*i+2] = best_val_2;
out.price = best_ent - best_val_1;
out.ent = best_ent;
out.val = best_val_0;
out.j = best_match;
return out;
}
int
auction_pass (const struct spcsr_t A, const double *b,
int *match_in,
double *price_in,
int n_unmatched,
int *unmatched_in,
const double mu,
const double feasbnd)
{
double least_val = HUGE_VAL;
int new_n_unmatched = 0;
int k;
int * restrict unmatched = unmatched_in;
int * restrict match = match_in;
double * restrict price = price_in;
for (k = 0; k < n_unmatched; ++k) {
const int i = unmatched[k];
struct auction_bid bid;
int old_i;
bid = auction_find_bid (A, b, i, price);
if (!isinf(bid.price) && bid.val < least_val) least_val = bid.val;
/*
if (bid.j == 182) {
fprintf (stderr, "bid: row %d for col %d, ent %g, price %g, val %g\n",
i, bid.j, bid.ent, bid.price, bid.val);
};
*/
assert (bid.price + mu >= price[bid.j]);
/*
if (1 == bid.j)
printf ("bidding for %d: (%d, %d), price %g\n", 1, i, bid.j, bid.price+mu);
*/
/* update matching */
old_i = match[bid.j];
assert(old_i != i);
if (old_i >= 0) {
/*
if (1 == old_i) printf ("unmatching (%d, %d) for (%d, %d)\n", old_i, bid.j, i, bid.j);
*/
unmatched[new_n_unmatched] = old_i;
++new_n_unmatched;
}
/*
if (1 == i) printf ("matching (%d, %d), price %g\n", i, bid.j, bid.price+mu);
*/
match[bid.j] = i;
price[bid.j] = bid.price + mu;
}
if (least_val < feasbnd) return -1;
return new_n_unmatched;
}
int
auction_pass_lifo (const struct spcsr_t A, const double *b,
int *match_in,
double *price_in,
int n_unmatched,
int *unmatched_in,
const double mu,
const double feasbnd)
{
double least_val = HUGE_VAL;
int new_n_unmatched = 0;
int k;
int * restrict unmatched = unmatched_in;
int * restrict match = match_in;
double * restrict price = price_in;
for (k = 0; k < n_unmatched; ++k) {
int i = unmatched[k];
struct auction_bid bid;
int old_i;
while (i >= 0) {
bid = auction_find_bid (A, b, i, price);
if (!isinf(bid.price) && bid.val < least_val) least_val = bid.val;
/*
if (bid.j == 182) {
fprintf (stderr, "bid: row %d for col %d, ent %g, price %g, val %g\n",
i, bid.j, bid.ent, bid.price, bid.val);
};
*/
assert (bid.price + mu >= price[bid.j]);
/*
if (1 == bid.j)
printf ("bidding for %d: (%d, %d), price %g\n", 1, i, bid.j, bid.price+mu);
*/
old_i = match[bid.j];
match[bid.j] = i;
price[bid.j] = bid.price + mu;
/* update matching */
assert(old_i != i);
i = old_i;
/*
if (1 == i) printf ("matching (%d, %d), price %g\n", i, bid.j, bid.price+mu);
*/
}
}
if (least_val < feasbnd) return -1;
return new_n_unmatched;
}
int
auction_pass_cache (const struct spcsr_t A, const double *b,
int *match_in,
double *price_in,
int n_unmatched,
int *unmatched_in,
int *jcache, double *ecache,
const double mu,
const double feasbnd)
{
double least_val = HUGE_VAL;
int new_n_unmatched = 0;
int k;
int * restrict unmatched = unmatched_in;
int * restrict match = match_in;
double * restrict price = price_in;
for (k = 0; k < n_unmatched; ++k) {
const int i = unmatched[k];
struct auction_bid bid;
int old_i;
bid = auction_find_bid_cache (A, b, i, price, jcache, ecache);
if (!isinf(bid.price) && bid.val < least_val) least_val = bid.val;
/*
if (bid.j == 182) {
fprintf (stderr, "bid: row %d for col %d, ent %g, price %g, val %g\n",
i, bid.j, bid.ent, bid.price, bid.val);
};
*/
assert (bid.price + mu >= price[bid.j]);
/*
if (1 == bid.j)
printf ("bidding for %d: (%d, %d), price %g\n", 1, i, bid.j, bid.price+mu);
*/
/* update matching */
old_i = match[bid.j];
assert(old_i != i);
if (old_i >= 0) {
/*
if (1 == old_i) printf ("unmatching (%d, %d) for (%d, %d)\n", old_i, bid.j, i, bid.j);
*/
unmatched[new_n_unmatched] = old_i;
++new_n_unmatched;
}
/*
if (1 == i) printf ("matching (%d, %d), price %g\n", i, bid.j, bid.price+mu);
*/
match[bid.j] = i;
price[bid.j] = bid.price + mu;
}
if (least_val < feasbnd) return -1;
return new_n_unmatched;
}
int
auction_pass_cl (const struct spcsr_t A, const double *b,
int *match_in,
double *price_in,
int n_unmatched,
int *unmatched_in,
int *n_changed_in,
int *changed_col_list_in,
int *col_changed_flag_in,
const double mu,
const double feasbnd)
{
double least_val = HUGE_VAL;
int new_n_unmatched = 0;
int k;
int * restrict unmatched = unmatched_in;
int * restrict match = match_in;
double * restrict price = price_in;
int n_changed = *n_changed_in;
int * restrict changed_col_list = changed_col_list_in;
int * restrict col_changed_flag = col_changed_flag_in;
for (k = 0; k < n_unmatched; ++k) {
const int i = unmatched[k];
struct auction_bid bid;
int old_i;
/*
if (i < 0) {
extern int rank;
fprintf (stderr, "%d: pass_cl n_un %d\n", rank, n_unmatched);
}
*/
assert(i >= 0);
bid = auction_find_bid (A, b, i, price);
if (!isinf(bid.price) && bid.val < least_val) least_val = bid.val;
#if 0
if (bid.j == 2050 || i == 680 || i == 2518 || i == 881) {
fprintf (stderr, "bid: row %d for col %d, ent %g, price %g, val %g\n",
i, bid.j, bid.ent, bid.price, bid.val);
};
#endif
/*
if (!(bid.price + mu >= price[bid.j])) {
extern int rank;
fprintf (stderr, "%d: %d: bid.price %g + mu %g < price[%d] %g\n",
rank, i, bid.price, mu, bid.j, price[bid.j]);
}
*/
assert (bid.price + mu >= price[bid.j]);
/*
if (1 == bid.j)
printf ("bidding for %d: (%d, %d), price %g\n", 1, i, bid.j, bid.price+mu);
*/
/* update matching */
old_i = match[bid.j];
assert(old_i != i);
if (old_i >= 0) {
/*
if (1 == old_i) printf ("unmatching (%d, %d) for (%d, %d)\n", old_i, bid.j, i, bid.j);
*/
unmatched[new_n_unmatched] = old_i;
++new_n_unmatched;
}
/*
if (1 == i) printf ("matching (%d, %d), price %g\n", i, bid.j, bid.price+mu);
*/
match[bid.j] = i;
price[bid.j] = bid.price + mu;
/* mark as changed */
if (col_changed_flag[bid.j] < 0) {
col_changed_flag[bid.j] = n_changed;
changed_col_list[n_changed] = bid.j;
++n_changed;
}
}
*n_changed_in = n_changed;
if (least_val < feasbnd) return -1;
return new_n_unmatched;
}
int lastiter;
inline int
auction_matchall (const struct spcsr_t A, const double *b,
int *match,
double *price,
int n_unmatched,
int *unmatched,
const double mu,
const double feasbnd)
{
int new_n_unmatched;
int cnt = 0;
while (n_unmatched > 0) {
/* printf ("... match[%d] = %d\n", 1, match[1]); */
#if 1
new_n_unmatched = auction_pass (A, b, match, price, n_unmatched, unmatched,
mu, feasbnd);
#elif 1
new_n_unmatched = auction_pass_lifo (A, b, match, price, n_unmatched, unmatched,
mu, feasbnd);
#else
{
#define N_CHUNK 10
int *unmatched_chunk = unmatched;
int n_unmatched_chunk, n_done = 0;
while (n_done < n_unmatched) {
unmatched_chunk = unmatched + n_done;
n_unmatched_chunk = N_CHUNK;
if (n_unmatched_chunk > n_unmatched - n_done)
n_unmatched_chunk = n_unmatched - n_done;
while (n_unmatched_chunk > 0)
n_unmatched_chunk = auction_pass (A, b, match, price,
n_unmatched_chunk, unmatched_chunk,
mu, feasbnd);
if (n_unmatched_chunk < 0) { n_unmatched = n_unmatched_chunk; break; }
n_done += N_CHUNK;
}
new_n_unmatched = 0;
}
#endif
n_unmatched = new_n_unmatched;
#if 0
if (new_n_unmatched == n_unmatched) ++cnt;
if (cnt > 3 && !lastiter) return new_n_unmatched;
#elif 0
if (new_n_unmatched == n_unmatched)
{
int j;
double mn = HUGE_VAL;
for (j = 0; j < A.nc; ++j) {
if (match[j] >= 0 && price[j] < mn) mn = price[j];
}
fprintf(stderr, "new_n_unmatched == n_unmatched == %d, mn %g\n", n_unmatched, mn);
if (mn != HUGE_VAL)
for (j = 0; j < A.nc; ++j)
if (match[j] < 0 && price[j] > mn) price[j] = mn;
}
#elif 0
if (new_n_unmatched == n_unmatched)
{
int k, j;
double mn = HUGE_VAL, e;
for (k = 0; k < A.nent; ++k) {
j = A.colind[k];
e = b[k];
if (match[j] < 0 && e < mn) mn = e;
}
fprintf(stderr, "new_n_unmatched == n_unmatched == %d, mn %g\n", n_unmatched, mn);
if (mn < HUGE_VAL)
for (k = 0; k < A.nc; ++k) {
if (match[k] < 0) price[k] += mn;
else price[k] -= mn;
}
}
#elif 0
if (new_n_unmatched == n_unmatched)
{
double mn = HUGE_VAL;
int k;
for (k = 0; k < new_n_unmatched; ++k) {
int i, k2;
i = unmatched[k];
for (k2 = A.rowoff[i]; k2 != A.rowoff[i+1]; ++k2) {
int j; double e;
j = A.colind[k2];
e = b[k2];
if (match[j] < 0 && e < mn) mn = e;
}
}
fprintf(stderr, "new_n_unmatched == n_unmatched == %d, mn %g\n", n_unmatched, mn);
if (mn < HUGE_VAL)
for (k = 0; k < A.nc; ++k) {
if (match[k] < 0) price[k] -= mn;
else price[k] += mn;
}
}
#elif 0
if (n_unmatched < 5) {
int k;
fprintf (stderr, "unmatched (%d):", n_unmatched);
for (k = 0; k < n_unmatched; ++k) {
fprintf (stderr, " %d", unmatched[k]);
}
fprintf (stderr, "\n");
}
#endif
}
return n_unmatched;
}
inline int
auction_matchall_adaptive (const struct spcsr_t A, const double *b,
int *match,
double *price,
int n_unmatched,
int *unmatched,
double *mu_in,
const double mu_top,
const double feasbnd)
{
int new_n_unmatched;
double mu = *mu_in;
int n_no_change = 0;
while (n_unmatched > 0) {
/* printf ("... match[%d] = %d\n", 1, match[1]); */
new_n_unmatched = auction_pass (A, b, match, price, n_unmatched, unmatched,
mu, feasbnd);
if (new_n_unmatched == n_unmatched) ++n_no_change;
if (n_no_change > 10) {
mu *= 2.0;
if (mu > mu_top) mu = mu_top;
n_no_change = 0;
}
n_unmatched = new_n_unmatched;
}
*mu_in = mu;
return n_unmatched;
}
inline int
auction_matchall_cache (const struct spcsr_t A, const double *b,
int *match,
double *price,
int n_unmatched,
int *unmatched,
int *jcache, double *ecache,
const double mu,
const double feasbnd)
{
int new_n_unmatched;
while (n_unmatched > 0) {
/* printf ("... match[%d] = %d\n", 1, match[1]); */
new_n_unmatched = auction_pass_cache (A, b, match, price,
n_unmatched, unmatched,
jcache, ecache,
mu, feasbnd);
n_unmatched = new_n_unmatched;
}
return n_unmatched;
}
int
auction_simple (const struct spcsr_t A, const double *b_in,
int *match, double *price,
const double mu_min)
{
int *unmatched = NULL;
int n_unmatched = 0;
double feasbnd;
if (A.nr > A.nc) return -1;
STATS(a_nstats = 0;);
memset (match, -1, A.nc * sizeof(int));
memset (price, 0, A.nc * sizeof(double));
unmatched = malloc (A.nr * sizeof(int));
for (n_unmatched = 0; n_unmatched < A.nr; ++n_unmatched)
unmatched[n_unmatched] = n_unmatched;
{
double mx = -HUGE_VAL, mn = HUGE_VAL;
int k;
const double * restrict b = b_in;
for (k = 0; k < A.nent; ++k) {
const double e = b[k];
if (e < mn) mn = e;
if (e > mx) mx = e;
}
feasbnd = -(2*A.nc-1) * (mx - mn) - (A.nc - 1) * mu_min;
}
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = HUGE_VAL;);
STATS(++a_nstats;);
lastiter = 1;
n_unmatched = auction_matchall (A, b_in, match, price,
n_unmatched, unmatched,
mu_min, feasbnd);
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = mu_min;);
STATS(auction_eval_gap(A, b_in, match, price,
&a_stats_gap[a_nstats]););
STATS(++a_nstats;);
free (unmatched);
return n_unmatched;
}
int
auction_simple_cache (const struct spcsr_t A, const double *b_in,
int *match, double *price,
const double mu_min)
{
int *unmatched = NULL;
int n_unmatched = 0;
double feasbnd;
int *jcache = NULL;
double *ecache = NULL;
if (A.nr > A.nc) return -1;
STATS(a_nstats = 0;);
memset (match, -1, A.nc * sizeof(int));
memset (price, 0, A.nc * sizeof(double));
unmatched = malloc (A.nr * sizeof(int));
for (n_unmatched = 0; n_unmatched < A.nr; ++n_unmatched)
unmatched[n_unmatched] = n_unmatched;
jcache = calloc(2 * A.nr, sizeof(int));
ecache = calloc(3 * A.nr, sizeof(double));
{
double mx = -HUGE_VAL, mn = HUGE_VAL;
int k;
const double * restrict b = b_in;
for (k = 0; k < A.nent; ++k) {
const double e = b[k];
if (e < mn) mn = e;
if (e > mx) mx = e;
}
feasbnd = -(2*A.nc-1) * (mx - mn) - (A.nc - 1) * mu_min;
}
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = HUGE_VAL;);
STATS(++a_nstats;);
n_unmatched = auction_matchall_cache (A, b_in, match, price,
n_unmatched, unmatched,
jcache, ecache,
mu_min, feasbnd);
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = mu_min;);
STATS(auction_eval_gap(A, b_in, match, price,
&a_stats_gap[a_nstats]););
STATS(++a_nstats;);
free (ecache);
free (jcache);
free (unmatched);
return n_unmatched;
}
int
auction_greedy_gather0 (const struct spcsr_t A, const double * restrict b,
int * restrict match,
double * restrict price,
int * restrict unmatched,
double * restrict bmn_out,
double * restrict bmx_out,
double * restrict gap_out,
double * restrict pmax_out)
{
int i, n_unmatched = 0;
double bmx = -HUGE_VAL, bmn = HUGE_VAL;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
double gap = 0.0;
/*
double primal = 0.0, dual = 0.0;
*/
double pmax = 0.0;
for (i = 0; i < A.nr; ++i) {
int k, kend = rowoff[i+1];
double best_e = -HUGE_VAL;
double largest_e = -HUGE_VAL;
int best_j = -1;
for (k = rowoff[i]; k < kend; ++k) {
const int j = colind[k];
const double e = b[k];
if (e < bmn) bmn = e;
if (e > bmx) bmx = e;
if (e > largest_e) largest_e = e;
if (e > best_e && match[j] < 0) {
best_e = e;
best_j = j;
}
}
/*dual += largest_e;*/ /* prices are zero */
if (best_j >= 0) {
/*primal += best_e;*/
gap += (largest_e - best_e);
price[best_j] = best_e;
if (best_e > pmax) pmax = best_e;
match[best_j] = i;
/* printf ("greedily matching (%d, %d)\n", i, best_j); */
} else {
gap += largest_e; /* an empty row will set gap to -inf */
/* printf ("greedily NOT matching %d [%d]\n", i, n_unmatched); */
unmatched[n_unmatched] = i;
++n_unmatched;
}
}
*bmx_out = bmx;
*bmn_out = bmn;
/*
*primal_out = primal;
*dual_out = dual;
*/
*gap_out = gap;
*pmax_out = pmax;
return n_unmatched;
}
int
auction_greedy_gather0_cl (const struct spcsr_t A, const double * restrict b,
int * restrict match,
double * restrict price,
int * restrict unmatched,
int * restrict n_changed_out,
int * restrict changed_col_list,
int * restrict col_changed_flag,
double * restrict bmn_out,
double * restrict bmx_out,
double * restrict gap_out,
double * restrict pmax_out)
{
int i, n_unmatched = 0;
double bmx = -HUGE_VAL, bmn = HUGE_VAL;
int n_changed = *n_changed_out;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
double gap = 0.0;
/*
double primal = 0.0, dual = 0.0;
*/
double pmax = 0.0;
for (i = 0; i < A.nr; ++i) {
int k, kend = rowoff[i+1];
double best_e = -HUGE_VAL;
double largest_e = -HUGE_VAL;
int best_j = -1;
for (k = rowoff[i]; k < kend; ++k) {
const int j = colind[k];
const double e = b[k];
if (e < bmn) bmn = e;
if (e > bmx) bmx = e;
if (e > largest_e) largest_e = e;
if (e > best_e && match[j] < 0) {
best_e = e;
best_j = j;
}
}
/*dual += largest_e;*/ /* prices are zero */
if (best_j >= 0) {
/*primal += best_e;*/
gap += (largest_e - best_e);
price[best_j] = best_e;
if (best_e > pmax) pmax = best_e;
match[best_j] = i;
/* printf ("greedily matching (%d, %d)\n", i, best_j); */
changed_col_list[n_changed] = best_j;
col_changed_flag[best_j] = n_changed;
++n_changed;
} else {
gap += largest_e; /* an empty row will set gap to -inf */
/* printf ("greedily NOT matching %d [%d]\n", i, n_unmatched); */
unmatched[n_unmatched] = i;
++n_unmatched;
}
}
*n_changed_out = n_changed;
*bmx_out = bmx;
*bmn_out = bmn;
/*
*primal_out = primal;
*dual_out = dual;
*/
*gap_out = gap;
*pmax_out = pmax;
return n_unmatched;
}
static int
greedy_unmatch (const struct spcsr_t A, const double * restrict b,
int * restrict match,
double * restrict price,
int n_unmatched,
int * restrict unmatched,
const double old_mu,
const double new_mu,
double * restrict pmax_out)
{
int j;
const int * restrict rowoff = A.rowoff;
const int * restrict colind = A.colind;
double pmax = 0.0;
for (j = 0; j < A.nc; ++j) {
const double p = price[j];
const int i = match[j];
double profit = -HUGE_VAL;
double match_e = -HUGE_VAL;
int k;
assert(i < A.nr);
if (i < 0) continue;
if (HUGE_VAL == p) continue;
if (p > pmax) pmax = p;
for (k = rowoff[i]; k != rowoff[i+1]; ++k) {
const int j2 = colind[k];
const double e = b[k];
double v;
v = e - p;
if (v > profit) profit = v;
if (j2 == j) {
match_e = e;
}
}
if ((p + profit) - match_e > new_mu) {
unmatched[n_unmatched] = i;
++n_unmatched;
match[j] = -1;
}
/*
else {
printf ("(%d,%d) price: %g profit: %g e: %g new_mu: %g\n diff %g\n",
i,j , p, profit, match_e, new_mu,
(p + profit) - match_e);
}
*/
#if 0
else {
price[j] += new_mu; /* place a bid... */
}
#endif
}
*pmax_out = pmax;
return n_unmatched;
}
static int
nomatch_gather0 (const struct spcsr_t A, const double * restrict b,
int * restrict unmatched,
double * restrict brange_out,
/*
double * restrict primal_out, double * restrict dual_out,
*/
double * restrict gap_out,
double * restrict pmax_out)
{
int i, n_unmatched = 0;
double bmx = -HUGE_VAL, bmn = HUGE_VAL;
const int * restrict rowoff = A.rowoff;
double brange, gap = 0.0;
/*double primal = 0.0, dual = 0.0;*/
for (i = 0; i < A.nr; ++i) {
int k, kend = rowoff[i+1];
double largest_e = -HUGE_VAL;
for (k = rowoff[i]; k < kend; ++k) {
const double e = b[k];
if (e < bmn) bmn = e;
if (e > bmx) bmx = e;
if (e > largest_e) largest_e = e;
}
gap += largest_e; /* prices are zero, no matching at all */
unmatched[n_unmatched] = i;
++n_unmatched;
}
brange = bmx - bmn;
*brange_out = brange;
/*
*primal_out = primal;
*dual_out = dual;
*/
*gap_out = gap;
*pmax_out = 0.0;
return n_unmatched;
}
int
auction_scaling (const struct spcsr_t A, const double *b_in,
int *match_in,
double *price_in,
const double mu_final, const int relgap)
{
int n_unmatched;
int *unmatched;
double mu, pmax, feasbnd;
double brange, bmn, bmx, gap;
/*double primal, dual;*/
const double * restrict b = b_in;
int * restrict match = match_in;
double * restrict price = price_in;
const double theta = 1.0/8.0;
#if !defined(NDEBUG)
int *seenr = NULL;
seenr = malloc (A.nr * sizeof(int));
memset (seenr, -1, A.nr * sizeof(int));
#endif
STATS(a_nstats = 0;);
n_unmatched = 0;
unmatched = malloc (A.nr * sizeof(int));
memset (match, -1, A.nc * sizeof(int));
memset (price, 0, A.nc * sizeof(double));
STATS(a_stats_n_unmatched[a_nstats] = A.nr;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = HUGE_VAL;);
STATS(++a_nstats;);
/* printf ("mu_final = %e\n", mu_final); */
#if 1
/* initial greedy matching while scanning largest benefit */
n_unmatched = auction_greedy_gather0 (A, b, match, price, unmatched,
&bmn, &bmx, &gap, &pmax);
brange = bmx - bmn;
mu = theta * gap / A.nc;
#else
n_unmatched = nomatch_gather0 (A, b, unmatched,
&brange, &gap, &pmax);
mu = brange / 16.0;
#endif
if (!n_unmatched) goto freemem;
if (gap == -HUGE_VAL) { n_unmatched = -1; goto freemem; }
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = gap;);
STATS(++a_nstats;);
if (mu < 8.0*mu_final) mu = mu_final;
feasbnd = -(2 * A.nc - 1) * brange - A.nc * mu_final;
for (;;) {
double new_mu;
if (!n_unmatched && mu <= mu_final)
break;
#if !defined(NDEBUG)
{ /* sanity checks */
int un;
for (un = 0; un < n_unmatched; ++un) {
const int i = unmatched[un];
int found = 0;
int j;
if (seenr[i] >= 0)
printf ("%d is duplicated in unmatched, %d %d\n", i, seenr[i], un);
else seenr[i] = un;
for (j = 0; j < A.nc; ++j) {
found = (match[j] == i);
if (found)
printf ("match (%d, %d) should have been cleared\n", i, j);
}
assert(!found);
}
for (un = 0; un < n_unmatched; ++un) seenr[unmatched[un]] = -1;
}
#endif
/* printf ("mu = %e, gap = %g\n", mu, dual - primal);*/
/* printf ("match[%d] = %d\n", 1, match[1]); */
#if 1
lastiter = (mu <= mu_final);
n_unmatched = auction_matchall (A, b, match_in, price_in,
n_unmatched, unmatched,
mu, feasbnd - pmax);
#elif 1
{
int n_no_change = 0;
while (n_unmatched > 0) {
int new_n_unmatched;
/*
fprintf (stderr, "%d: n_unmatched before pass_cl: %d\n", rank,
n_unmatched);
*/
new_n_unmatched = auction_pass (A, b, match, price,
n_unmatched, unmatched,
mu, -HUGE_VAL);
if (new_n_unmatched == n_unmatched) {
++n_no_change;
if (n_no_change > 5) {
n_unmatched = new_n_unmatched;
break;
}
}
else
n_no_change = 0;
n_unmatched = new_n_unmatched;
/*
fprintf (stderr, "%d: n_unmatched after pass_cl: %d\n", rank,
n_unmatched);
*/
}
}
#else
{ /* doesn't seem to buy _anything_ */
double mu_a = mu_final;
n_unmatched = auction_matchall_adaptive (A, b, match_in, price_in,
n_unmatched, unmatched,
&mu_a, mu, feasbnd - pmax);
fprintf(stderr, "adaptive: %g, top %g, bot %g\n", mu_a, mu,
mu_final);
mu = mu_a;
}
#endif
if (n_unmatched < 0 || mu <= mu_final)
break;
new_mu = mu;
{
int j;
pmax = -HUGE_VAL;
if (!relgap) {
auction_eval_gap (A, b, match, price, &gap);
#if 1
if (gap < A.nc * mu_final) break; /* before clearing match */
#else
if (gap < mu_final) break; /* before clearing match */
#endif
}
else {
double pval, dval;
auction_eval_primal_mdual (A, b, match, price, &pval, &dval);
gap = dval - pval;
if (gap / pval < mu_final) break;
}
new_mu = theta * (gap) / A.nc;
if (new_mu < 8.0*mu_final) new_mu = mu_final;
if (new_mu > theta * mu) new_mu = theta * mu;
#if 1
for (j = 0; j < A.nc; ++j) {
const int i = match[j];
const double p = price[j];
if (p != HUGE_VAL) {
/* if (1 == i) printf ("clearing (%d, %d)\n", i, j); */
match[j] = -1;
unmatched[n_unmatched] = i;
++n_unmatched;
if (p > pmax) pmax = p;
}
/*
else
if (1 == i) printf ("link (%d, %d) is permanent\n", i, j);
*/
}
#elif 1
n_unmatched = greedy_unmatch (A, b, match, price,
n_unmatched, unmatched,
mu, mu_final, /*0 or new_mu,*/
&pmax);
if (!n_unmatched) break;
#else
for (;;)
{
int new_n_unmatched;
new_n_unmatched = greedy_unmatch (A, b, match, price,
n_unmatched, unmatched,
mu, new_mu, /*0 or new_mu,*/
&pmax);
/*fprintf (stderr, "after gr unmatch: %d\n", n_unmatched);*/
if (new_n_unmatched == n_unmatched) break;
n_unmatched = new_n_unmatched;
new_mu /= 4.0;
}
/*if (!n_unmatched) break;*/
#endif
/*
fprintf (stderr, "after unmatch: %d\n", n_unmatched);
*/
}
/*
printf( "... mu = %g\n... n_unmatched = %d\n", mu, n_unmatched);
printf( "... gap = %g\n", gap);
*/
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = mu;);
STATS(a_stats_gap[a_nstats] = gap;);
STATS(++a_nstats;);
mu = new_mu;
}
#if !defined(NOSTATS)
a_stats_n_unmatched[a_nstats] = n_unmatched;
a_stats_mu[a_nstats] = mu;
auction_eval_gap (A, b, match, price,
&a_stats_gap[a_nstats]);
++a_nstats;
#endif
freemem:
#if !defined(NDEBUG)
free (seenr);
#endif
free (unmatched);
return n_unmatched;
}
int
auction_scaling_cache (const struct spcsr_t A, const double *b_in,
int *match_in,
double *price_in,
const double mu_final)
{
int n_unmatched;
int *unmatched;
double mu, pmax, feasbnd;
double brange, gap;
/*double primal, dual;*/
const double * restrict b = b_in;
int * restrict match = match_in;
double * restrict price = price_in;
int *jcache = NULL;
double *ecache = NULL;
const double theta = 1.0/8.0;
#if !defined(NDEBUG)
int *seenr = NULL;
seenr = malloc (A.nr * sizeof(int));
memset (seenr, -1, A.nr * sizeof(int));
#endif
STATS(a_nstats = 0;);
jcache = calloc(2 * A.nr, sizeof(int));
ecache = calloc(3 * A.nr, sizeof(double));
n_unmatched = 0;
unmatched = malloc (A.nr * sizeof(int));
memset (match, -1, A.nc * sizeof(int));
memset (price, 0, A.nc * sizeof(double));
STATS(a_stats_n_unmatched[a_nstats] = A.nr;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = HUGE_VAL;);
STATS(++a_nstats;);
/* printf ("mu_final = %e\n", mu_final); */
#if 1
/* initial greedy matching while scanning largest benefit */
{ double bmn, bmx;
n_unmatched = auction_greedy_gather0 (A, b, match, price, unmatched,
&bmn, &bmx, &gap, &pmax);
brange = bmx - bmn;
}
mu = theta * gap / A.nc;
#else
n_unmatched = nomatch_gather0 (A, b, unmatched,
&brange, &gap, &pmax);
mu = brange / 16.0;
#endif
if (!n_unmatched) goto freemem;
if (gap == -HUGE_VAL) { n_unmatched = -1; goto freemem; }
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = HUGE_VAL;);
STATS(a_stats_gap[a_nstats] = gap;);
STATS(++a_nstats;);
if (mu < 8.0*mu_final) mu = mu_final;
feasbnd = -(2 * A.nc - 1) * brange - A.nc * mu_final;
for (;;) {
double new_mu;
#if !defined(NDEBUG)
{ /* sanity checks */
int un;
for (un = 0; un < n_unmatched; ++un) {
const int i = unmatched[un];
int found = 0;
int j;
if (seenr[i] >= 0)
printf ("%d is duplicated in unmatched, %d %d\n", i, seenr[i], un);
else seenr[i] = un;
for (j = 0; j < A.nc; ++j) {
found = (match[j] == i);
if (found)
printf ("match (%d, %d) should have been cleared\n", i, j);
}
assert(!found);
}
for (un = 0; un < n_unmatched; ++un) seenr[unmatched[un]] = -1;
}
#endif
/* printf ("mu = %e, gap = %g\n", mu, dual - primal);*/
/* printf ("match[%d] = %d\n", 1, match[1]); */
n_unmatched = auction_matchall_cache (A, b, match_in, price_in,
n_unmatched, unmatched,
jcache, ecache,
mu, feasbnd - pmax);
if (n_unmatched || mu <= mu_final)
break;
{
int j;
pmax = -HUGE_VAL;
auction_eval_gap (A, b, match, price, &gap);
#if 1
if (gap < A.nc * mu_final) break; /* before clearing match */
#else
if (gap < mu_final) break; /* before clearing match */
#endif
new_mu = theta * (gap) / A.nc;
if (new_mu < 8.0*mu_final) new_mu = mu_final;
if (new_mu > theta * mu) new_mu = theta * mu;
#if 1
for (j = 0; j < A.nc; ++j) {
const int i = match[j];
const double p = price[j];
if (p != HUGE_VAL) {
/* if (1 == i) printf ("clearing (%d, %d)\n", i, j); */
match[j] = -1;
unmatched[n_unmatched] = i;
++n_unmatched;
if (p > pmax) pmax = p;
}
/*
else
if (1 == i) printf ("link (%d, %d) is permanent\n", i, j);
*/
}
#else
n_unmatched = greedy_unmatch (A, b, match, price, unmatched,
mu, mu_final,
&pmax);
if (!n_unmatched) break;
#endif
}
/*
printf( "... mu = %g\n... n_unmatched = %d\n", mu, n_unmatched);
printf( "... gap = %g\n", gap);
*/
STATS(a_stats_n_unmatched[a_nstats] = n_unmatched;);
STATS(a_stats_mu[a_nstats] = mu;);
STATS(a_stats_gap[a_nstats] = gap;);
STATS(++a_nstats;);
mu = new_mu;
}
#if !defined(NOSTATS)
a_stats_n_unmatched[a_nstats] = n_unmatched;
a_stats_mu[a_nstats] = mu;
auction_eval_gap (A, b, match, price,
&a_stats_gap[a_nstats]);
++a_nstats;
#endif
freemem:
#if !defined(NDEBUG)
free (seenr);
#endif
free (ecache);
free (jcache);
free (unmatched);
return n_unmatched;
}
void
a_printstats (void)
{
#if !defined(NOSTATS)
int k;
printf("n_unmatched:");
for (k = 0; k < a_nstats; ++k) {
printf (" %d", a_stats_n_unmatched[k]);
}
printf ("\n");
printf("mu:");
for (k = 0; k < a_nstats; ++k) {
printf (" %g", a_stats_mu[k]);
}
printf ("\n");
printf("gap:");
for (k = 0; k < a_nstats; ++k) {
printf (" %g", a_stats_gap[k]);
}
printf ("\n");
#endif
}