/* * Set up environment and global constraints (dir-edge constraints, containment constraints * etc). * * diredges: 0=no dir edge constraints * 1=one separation constraint for each edge (in acyclic subgraph) * 2=DiG-CoLa level constraints */ CMajEnvVPSC *initCMajVPSC(int n, float *packedMat, vtx_data * graph, ipsep_options * opt, int diredges) { int i, j; /* nv is the number of real nodes */ int nConCs; /* fprintf(stderr,"Entered initCMajVPSC\n"); */ CMajEnvVPSC *e = GNEW(CMajEnvVPSC); e->A = NULL; e->packedMat = packedMat; /* if we have clusters then we'll need two constraints for each var in * a cluster */ e->nldv = 2 * opt->clusters->nclusters; e->nv = n - e->nldv; e->ndv = 0; e->gcs = NULL; e->vs = N_GNEW(n, Variable *); for (i = 0; i < n; i++) { e->vs[i] = newVariable(i, 1.0, 1.0); } e->gm = 0; if (diredges == 1) { if (Verbose) fprintf(stderr, " generate edge constraints...\n"); for (i = 0; i < e->nv; i++) { for (j = 1; j < graph[i].nedges; j++) { /* fprintf(stderr,"edist=%f\n",graph[i].edists[j]); */ if (graph[i].edists[j] > 0.01) { e->gm++; } } } e->gcs = newConstraints(e->gm); e->gm = 0; for (i = 0; i < e->nv; i++) { for (j = 1; j < graph[i].nedges; j++) { int u = i, v = graph[i].edges[j]; if (graph[i].edists[j] > 0) { e->gcs[e->gm++] = newConstraint(e->vs[u], e->vs[v], opt->edge_gap); } } } } else if (diredges == 2) { int *ordering = NULL, *ls = NULL, cvar; double halfgap; DigColaLevel *levels; Variable **vs = e->vs; /* e->ndv is the number of dummy variables required, one for each boundary */ if (compute_hierarchy(graph, e->nv, 1e-2, 1e-1, NULL, &ordering, &ls, &e->ndv)) return NULL; levels = assign_digcola_levels(ordering, e->nv, ls, e->ndv); if (Verbose) fprintf(stderr, "Found %d DiG-CoLa boundaries\n", e->ndv); e->gm = get_num_digcola_constraints(levels, e->ndv + 1) + e->ndv - 1; e->gcs = newConstraints(e->gm); e->gm = 0; e->vs = N_GNEW(n + e->ndv, Variable *); for (i = 0; i < n; i++) { e->vs[i] = vs[i]; } free(vs); /* create dummy vars */ for (i = 0; i < e->ndv; i++) { /* dummy vars should have 0 weight */ cvar = n + i; e->vs[cvar] = newVariable(cvar, 1.0, 0.000001); } halfgap = opt->edge_gap; for (i = 0; i < e->ndv; i++) { cvar = n + i; /* outgoing constraints for each var in level below boundary */ for (j = 0; j < levels[i].num_nodes; j++) { e->gcs[e->gm++] = newConstraint(e->vs[levels[i].nodes[j]], e->vs[cvar], halfgap); } /* incoming constraints for each var in level above boundary */ for (j = 0; j < levels[i + 1].num_nodes; j++) { e->gcs[e->gm++] = newConstraint(e->vs[cvar], e->vs[levels[i + 1].nodes[j]], halfgap); } } /* constraints between adjacent boundary dummy vars */ for (i = 0; i < e->ndv - 1; i++) { e->gcs[e->gm++] = newConstraint(e->vs[n + i], e->vs[n + i + 1], 0); } }
/********************** lap: the upper RHS of the symmetric graph laplacian matrix which will be transformed to the hessian of the non-linear part of the optimisation function n: number of nodes (length of coords array) ordering: array containing sequences of nodes for each level, ie, ordering[levels[i]] is first node of (i+1)th level level_indexes: array of starting node for each level in ordering ie, levels[i] is index to first node of (i+1)th level also, levels[0] is number of nodes in first level and, levels[i]-levels[i-1] is number of nodes in ith level and, n - levels[num_divisions-1] is number of nodes in last level num_divisions: number of divisions between levels, ie number of levels - 1 separation: the minimum separation between nodes on different levels ***********************/ MosekEnv *mosek_init_hier(float *lap, int n, int *ordering, int *level_indexes, int num_divisions, float separation) { int count = 0; int i, j, num_levels = num_divisions + 1; int num_constraints; MosekEnv *mskEnv = GNEW(MosekEnv); DigColaLevel *levels; int nonzero_lapsize = (n * (n - 1)) / 2; /* vars for nodes (except x0) + dummy nodes between levels * x0 is fixed at 0, and therefore is not included in the opt problem * add 2 more vars for top and bottom constraints */ mskEnv->num_variables = n + num_divisions + 1; logfile = fopen("quad_solve_log", "w"); levels = assign_digcola_levels(ordering, n, level_indexes, num_divisions); #ifdef DUMP_CONSTRAINTS print_digcola_levels(logfile, levels, num_levels); #endif /* nonlinear coefficients matrix of objective function */ /* int lapsize=mskEnv->num_variables+(mskEnv->num_variables*(mskEnv->num_variables-1))/2; */ mskEnv->qval = N_GNEW(nonzero_lapsize, double); mskEnv->qsubi = N_GNEW(nonzero_lapsize, int); mskEnv->qsubj = N_GNEW(nonzero_lapsize, int); /* solution vector */ mskEnv->xx = N_GNEW(mskEnv->num_variables, double); /* constraint matrix */ separation /= 2.0; /* separation between each node and it's adjacent constraint */ num_constraints = get_num_digcola_constraints(levels, num_levels) + num_divisions + 1; /* constraints of the form x_i - x_j >= sep so 2 non-zero entries per constraint in LHS matrix * except x_0 (fixed at 0) constraints which have 1 nz val each. */ #ifdef EQUAL_WIDTH_LEVELS num_constraints += num_divisions; #endif /* pointer to beginning of nonzero sequence in a column */ for (i = 0; i < n - 1; i++) { for (j = i; j < n - 1; j++) { mskEnv->qval[count] = -2 * lap[count + n]; assert(mskEnv->qval[count] != 0); mskEnv->qsubi[count] = j; mskEnv->qsubj[count] = i; count++; } } #ifdef DUMP_CONSTRAINTS fprintf(logfile, "Q=["); int lapcntr = n; for (i = 0; i < mskEnv->num_variables; i++) { if (i != 0) fprintf(logfile, ";"); for (j = 0; j < mskEnv->num_variables; j++) { if (j < i || i >= n - 1 || j >= n - 1) { fprintf(logfile, "0 "); } else { fprintf(logfile, "%f ", -2 * lap[lapcntr++]); } } } fprintf(logfile, "]\nQ=Q-diag(diag(Q))+Q'\n"); #endif fprintf(logfile, "\n"); /* Make the mosek environment. */ mskEnv->r = MSK_makeenv(&mskEnv->env, NULL, NULL, NULL, NULL); /* Check whether the return code is ok. */ if (mskEnv->r == MSK_RES_OK) { /* Directs the log stream to the user * specified procedure 'printstr'. */ MSK_linkfunctoenvstream(mskEnv->env, MSK_STREAM_LOG, NULL, printstr); } /* Initialize the environment. */ mskEnv->r = MSK_initenv(mskEnv->env); if (mskEnv->r == MSK_RES_OK) { /* Make the optimization task. */ mskEnv->r = MSK_maketask(mskEnv->env, num_constraints, mskEnv->num_variables, &mskEnv->task); if (mskEnv->r == MSK_RES_OK) { int c_ind = 0; int c_var = n - 1; mskEnv->r = MSK_linkfunctotaskstream(mskEnv->task, MSK_STREAM_LOG, NULL, printstr); /* Resize the task. */ if (mskEnv->r == MSK_RES_OK) mskEnv->r = MSK_resizetask(mskEnv->task, num_constraints, mskEnv->num_variables, 0, /* no cones!! */ /* each constraint applies to 2 vars */ 2 * num_constraints + num_divisions, nonzero_lapsize); /* Append the constraints. */ if (mskEnv->r == MSK_RES_OK) mskEnv->r = MSK_append(mskEnv->task, 1, num_constraints); /* Append the variables. */ if (mskEnv->r == MSK_RES_OK) mskEnv->r = MSK_append(mskEnv->task, 0, mskEnv->num_variables); /* Put variable bounds. */ for (j = 0; j < mskEnv->num_variables && mskEnv->r == MSK_RES_OK; ++j) mskEnv->r = MSK_putbound(mskEnv->task, 0, j, MSK_BK_RA, -MSK_INFINITY, MSK_INFINITY); for (j = 0; j < levels[0].num_nodes && mskEnv->r == MSK_RES_OK; j++) { int node = levels[0].nodes[j] - 1; if (node >= 0) { INIT_sub_val(c_var,node); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 2, subi, vali); } else { /* constraint for y0 (fixed at 0) */ mskEnv->r = MSK_putaij(mskEnv->task, c_ind, c_var, 1.0); } mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_LO, separation, MSK_INFINITY); c_ind++; } for (i = 0; i < num_divisions && mskEnv->r == MSK_RES_OK; i++) { c_var = n + i; for (j = 0; j < levels[i].num_nodes && mskEnv->r == MSK_RES_OK; j++) { /* create separation constraint a>=b+separation */ int node = levels[i].nodes[j] - 1; if (node >= 0) { /* no constraint for fixed node */ INIT_sub_val(node,c_var); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 2, subi, vali); } else { /* constraint for y0 (fixed at 0) */ mskEnv->r = MSK_putaij(mskEnv->task, c_ind, c_var, -1.0); } mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_LO, separation, MSK_INFINITY); c_ind++; } for (j = 0; j < levels[i + 1].num_nodes && mskEnv->r == MSK_RES_OK; j++) { int node = levels[i + 1].nodes[j] - 1; if (node >= 0) { INIT_sub_val(c_var,node); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 2, subi, vali); } else { /* constraint for y0 (fixed at 0) */ mskEnv->r = MSK_putaij(mskEnv->task, c_ind, c_var, 1.0); } mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_LO, separation, MSK_INFINITY); c_ind++; } } c_var = n + i; for (j = 0; j < levels[i].num_nodes && mskEnv->r == MSK_RES_OK; j++) { /* create separation constraint a>=b+separation */ int node = levels[i].nodes[j] - 1; if (node >= 0) { /* no constraint for fixed node */ INIT_sub_val(node,c_var); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 2, subi, vali); } else { /* constraint for y0 (fixed at 0) */ mskEnv->r = MSK_putaij(mskEnv->task, c_ind, c_var, -1.0); } mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_LO, separation, MSK_INFINITY); c_ind++; } /* create constraints preserving the order of dummy vars */ for (i = 0; i < num_divisions + 1 && mskEnv->r == MSK_RES_OK; i++) { int c_var = n - 1 + i, c_var2 = c_var + 1; INIT_sub_val(c_var,c_var2); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 2, subi, vali); mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_LO, 0, MSK_INFINITY); c_ind++; } #ifdef EQUAL_WIDTH_LEVELS for (i = 1; i < num_divisions + 1 && mskEnv->r == MSK_RES_OK; i++) { int c_var = n - 1 + i, c_var_lo = c_var - 1, c_var_hi = c_var + 1; INIT_sub_val3(c_var_lo, c_var, c_var_h); mskEnv->r = MSK_putavec(mskEnv->task, 1, c_ind, 3, subi, vali); mskEnv->r = MSK_putbound(mskEnv->task, 1, c_ind, MSK_BK_FX, 0, 0); c_ind++; } #endif assert(c_ind == num_constraints); #ifdef DUMP_CONSTRAINTS fprintf(logfile, "A=["); for (i = 0; i < num_constraints; i++) { if (i != 0) fprintf(logfile, ";"); for (j = 0; j < mskEnv->num_variables; j++) { double aij; MSK_getaij(mskEnv->task, i, j, &aij); fprintf(logfile, "%f ", aij); } } fprintf(logfile, "]\n"); fprintf(logfile, "b=["); for (i = 0; i < num_constraints; i++) { fprintf(logfile, "%f ", separation); } fprintf(logfile, "]\n"); #endif if (mskEnv->r == MSK_RES_OK) { /* * The lower triangular part of the Q * matrix in the objective is specified. */ mskEnv->r = MSK_putqobj(mskEnv->task, nonzero_lapsize, mskEnv->qsubi, mskEnv->qsubj, mskEnv->qval); } } } delete_digcola_levels(levels, num_levels); return mskEnv; }