3 /* Copyright (c) 2007 Arnaud Legrand, Pedro Velho. All rights reserved. */
5 /* This program is free software; you can redistribute it and/or modify it
6 * under the terms of the license (GNU LGPL) which comes with this package. */
10 #include "xbt/sysdep.h"
11 #include "xbt/mallocator.h"
12 #include "maxmin_private.h"
20 * SDP specific variables.
22 #include <declarations.h>
25 static void create_cross_link(struct constraintmatrix *myconstraints, int k);
27 static void addentry(struct constraintmatrix *constraints,
28 struct blockmatrix *, int matno, int blkno, int indexi, int indexj, double ent, int blocksize);
31 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(surf_sdp, surf,
32 "Logging specific to SURF (sdp)");
34 ########################################################################
35 ######################## Simple Proportionnal fairness #################
36 ########################################################################
37 ####### Original problem ##########
46 # We assume in the following that n=2^K
48 ####### Standard SDP form #########
50 # A SDP can be written in the following standard form:
52 # (P) min c1*x1+c2*x2+...+cm*xm
53 # st F1*x1+F2*x2+...+Fm*xm-F0=X
56 # Where F1, F2, ..., Fm are symetric matrixes of size n by n. The
57 # constraint X>0 means that X must be positive semidefinite positive.
59 ########## Equivalent SDP #########
63 # y(k,i) for k in [0,K] and i in [1,2^k]
68 # y(2,1) y(2,2) y(2,3) y(2,4)
69 # y(3,1) y(3,2) y(3,3) y(3,4) y(3,5) y(3,6) y(3,7) y(3,8)
70 # -------------------------------------------------------
71 # x_1 x_2 x_3 x_4 x_5 x_6 x_7 x_8
74 # Structure Constraints:
76 # forall k in [0,K-1], and i in [1,2^k]
77 # [ y(k, 2i-1) y(k-1, i) ]
78 # [ y(k-1, i ) y(k, 2i) ]
82 # Capacity Constraints:
84 # -(A.y(K,*))_l >= - b_l
86 # Positivity Constraints:
87 # forall k in [0,K], and i in [1,2^k]
93 //typedef struct lmm_system {
95 // s_xbt_swag_t variable_set; /* a list of lmm_variable_t */
96 // s_xbt_swag_t constraint_set; /* a list of lmm_constraint_t */
97 // s_xbt_swag_t active_constraint_set; /* a list of lmm_constraint_t */
98 // s_xbt_swag_t saturated_variable_set; /* a list of lmm_variable_t */
99 // s_xbt_swag_t saturated_constraint_set; /* a list of lmm_constraint_t_t */
100 // xbt_mallocator_t variable_mallocator;
103 #define get_y(a,b) (pow(2,a)+b-1)
105 void sdp_solve(lmm_system_t sys)
107 lmm_variable_t var = NULL;
110 * SDP mapping variables.
117 int nb_cnsts_capacity=0;
118 int nb_cnsts_struct=0;
119 int nb_cnsts_positivy=0;
123 FILE *sdpout = fopen("SDPA-printf.tmp","w");
125 double *tempdiag = NULL;
132 * CSDP library specific variables.
134 struct blockmatrix C;
135 struct blockmatrix X,Z;
139 struct constraintmatrix *constraints;
142 * Classic maxmin variables.
144 lmm_constraint_t cnst = NULL;
145 lmm_element_t elem = NULL;
146 xbt_swag_t cnst_list = NULL;
147 xbt_swag_t var_list = NULL;
148 xbt_swag_t elem_list = NULL;
150 if ( !(sys->modified))
156 * Initialize the var list variable with only the active variables.
157 * Associate an index in the swag variables.
160 var_list = &(sys->variable_set);
161 DEBUG1("Variable set : %d", xbt_swag_size(var_list));
162 xbt_swag_foreach(var, var_list) {
164 if(var->weight) var->index = i++;
166 cnst_list=&(sys->saturated_constraint_set);
169 * Those fields are the top level description of the platform furnished in the xml file.
172 links = xbt_swag_size(&(sys->active_constraint_set));
175 * This number is found based on the tree structure explained on top.
177 K = (int)log((double)flows)/log(2.0);
180 * The number of variables in the SDP style.
182 nb_var = get_y(K, pow(2,K));
185 * Find the size of each group of constraints.
187 nb_cnsts_capacity = links + ((int)pow(2,K)) - flows;
188 nb_cnsts_struct = (int)pow(2,K) - 1;
189 nb_cnsts_positivy = (int)pow(2,K);
192 * The total number of constraints.
194 nb_cnsts = nb_cnsts_capacity + nb_cnsts_struct + nb_cnsts_positivy;
198 * Keep track of which blocks have off diagonal entries.
200 isdiag=(int *)calloc((nb_cnsts+1), sizeof(int));
201 for (i=1; i<=nb_cnsts; i++)
204 C.nblocks = nb_cnsts;
205 C.blocks = (struct blockrec *) calloc((C.nblocks+1),sizeof(struct blockrec));
206 constraints = (struct constraintmatrix *)calloc((nb_var+1),sizeof(struct constraintmatrix));
208 for(i = 1; i <= nb_var; i++){
209 constraints[i].blocks=NULL;
212 a = (double *)calloc(nb_var+1, sizeof(double));
221 * For each constraint block do.
223 for(i = 1; i <= nb_cnsts; i++){
226 * Structured blocks are size 2 and all others are size 1.
228 if(i <= nb_cnsts_struct){
230 fprintf(sdpout,"2 ");
233 fprintf(sdpout,"1 ");
237 * All blocks are matrices.
239 C.blocks[block_num].blockcategory = MATRIX;
240 C.blocks[block_num].blocksize = block_size;
241 C.blocks[block_num].data.mat = (double *)calloc(block_size * block_size, sizeof(double));
246 fprintf(sdpout,"\n");
250 * Creates de objective function array.
252 a = (double *)calloc((nb_var+1), sizeof(double));
254 for(i = 1; i <= nb_var; i++){
256 fprintf(sdpout,"-1 ");
259 fprintf(sdpout,"0 ");
263 fprintf(sdpout,"\n");
267 * Structure contraint blocks.
271 for(k = 1; k <= K; k++){
272 for(i = 1; i <= pow(2,k-1); i++){
273 matno=get_y(k,2*i-1);
274 fprintf(sdpout,"%d %d 1 1 1\n", matno , block_num);
275 addentry(constraints, &C, matno, block_num, 1, 1, 1.0, C.blocks[block_num].blocksize);
278 fprintf(sdpout,"%d %d 2 2 1\n", matno , block_num);
279 addentry(constraints, &C, matno, block_num, 2, 2, 1.0, C.blocks[block_num].blocksize);
282 fprintf(sdpout,"%d %d 1 2 1\n", matno , block_num);
283 addentry(constraints, &C, matno, block_num, 1, 2, 1.0, C.blocks[block_num].blocksize);
286 fprintf(sdpout,"%d %d 2 1 1\n", matno , block_num);
287 addentry(constraints, &C, matno, block_num, 2, 1, 1.0, C.blocks[block_num].blocksize);
289 isdiag[block_num] = 0;
296 * Capacity constraint block.
298 xbt_swag_foreach(cnst, cnst_list) {
300 fprintf(sdpout,"0 %d 1 1 %d\n", block_num, (int) - (cnst->bound));
301 addentry(constraints, &C, 0, block_num, 1, 1, - (cnst->bound) , C.blocks[block_num].blocksize);
303 elem_list = &(cnst->element_set);
304 xbt_swag_foreach(elem, elem_list) {
305 if(elem->variable->weight <=0) break;
306 fprintf(sdpout,"%d %d 1 1 %d\n", elem->variable->index, block_num, (int) - (elem->variable->value));
307 addentry(constraints, &C, elem->variable->index, block_num, 1, 1, - (elem->value), C.blocks[block_num].blocksize);
312 //Positivy constraint blocks
313 for(i = 1; i <= pow(2,K); i++){
315 fprintf(sdpout,"%d %d 1 1 1\n", matno, block_num);
316 addentry(constraints, &C, matno, block_num, 1, 1, 1.0, C.blocks[block_num].blocksize);
322 * At this point, we'll stop to recognize whether any of the blocks
323 * are "hidden LP blocks" and correct the block type if needed.
325 for (i=1; i<=nb_cnsts; i++){
326 if ((C.blocks[i].blockcategory != DIAG) &&
327 (isdiag[i]==1) && (C.blocks[i].blocksize > 1)){
329 * We have a hidden diagonal block!
332 printf("Block %d is actually diagonal.\n",i);
334 blocksz=C.blocks[i].blocksize;
335 tempdiag=(double *)calloc((blocksz+1), sizeof(double));
336 for (j=1; j<=blocksz; j++)
337 tempdiag[j]=C.blocks[i].data.mat[ijtok(j,j,blocksz)];
338 free(C.blocks[i].data.mat);
339 C.blocks[i].data.vec=tempdiag;
340 C.blocks[i].blockcategory=DIAG;
346 * Next, setup issparse and NULL out all nextbyblock pointers.
348 struct sparseblock *p=NULL;
349 for (i=1; i<=k; i++) {
350 p=constraints[i].blocks;
353 * First, set issparse.
355 if (((p->numentries) > 0.25*(p->blocksize)) && ((p->numentries) > 15)){
361 if (C.blocks[p->blocknum].blockcategory == DIAG)
365 * Setup the cross links.
375 * Create cross link reference.
377 create_cross_link(constraints, nb_var);
381 * Debuging print problem in SDPA format.
383 printf("Printing SDPA...\n");
384 if(XBT_LOG_ISENABLED(surf_sdp, xbt_log_priority_debug)) {
385 char *tmp=strdup("SDPA.tmp");
386 write_prob(tmp, nb_cnsts, nb_var, C, a, constraints);
387 //int write_prob(char *fname, int n, int k, struct blockmatrix C, double *a, struct constraintmatrix *constraints);
392 * Initialize parameters.
394 printf("Initializing solution...\n");
395 initsoln(nb_cnsts, nb_var, C, a, constraints, &X, &y, &Z);
402 printf("Calling the solver...\n");
403 int ret = easy_sdp(nb_cnsts, nb_var, C, a, constraints, 0.0, &X, &y, &Z, &pobj, &dobj);
407 case 1: printf("SUCCESS The problem is primal infeasible\n");
410 case 2: printf("SUCCESS The problem is dual infeasible\n");
413 case 3: printf("Partial SUCCESS A solution has been found, but full accuracy was not achieved. One or more of primal infeasibility, dual infeasibility, or relative duality gap are larger than their tolerances, but by a factor of less than 1000.\n");
416 case 4: printf("Failure. Maximum number of iterations reached.");
419 case 5: printf("Failure. Stuck at edge of primal feasibility.");
425 * Write out the solution if necessary.
427 //printf("Writting simple dsp...\n");
428 //write_sol("output.sol", n, k, X, y, Z);
433 //free_prob(n, k, C, a, constraints, X, y, Z);
439 if(XBT_LOG_ISENABLED(surf_sdp, xbt_log_priority_debug)) {
447 * Create the cross_link reference in order to have a byblock list.
449 void create_cross_link(struct constraintmatrix *myconstraints, int k){
453 struct sparseblock *p;
454 struct sparseblock *q;
456 struct sparseblock *prev;
462 for (i=1; i<=k; i++){
463 p=myconstraints[i].blocks;
465 if (p->nextbyblock == NULL){
469 * link in the remaining blocks.
471 for (j=i+1; j<=k; j++){
472 q=myconstraints[j].blocks;
475 if (q->blocknum == p->blocknum){
476 if (p->nextbyblock == NULL){
500 void addentry(struct constraintmatrix *constraints,
501 struct blockmatrix *C,
509 struct sparseblock *p;
510 struct sparseblock *p_sav;
512 p=constraints[matno].blocks;
517 * We haven't yet allocated any blocks.
519 p=(struct sparseblock *)calloc(1, sizeof(struct sparseblock));
521 //two entries because this library ignores indices starting in zerox
522 p->constraintnum=matno;
527 p->entries=calloc(p->numentries+1, sizeof(double));
528 p->iindices=calloc(p->numentries+1, sizeof(int));
529 p->jindices=calloc(p->numentries+1, sizeof(int));
531 p->entries[p->numentries]=ent;
532 p->iindices[p->numentries]=indexi;
533 p->jindices[p->numentries]=indexj;
535 p->blocksize=blocksize;
537 constraints[matno].blocks=p;
540 * We have some existing blocks. See whether this block is already
544 if (p->blocknum == blkno){
546 * Found the right block.
548 p->constraintnum=matno;
550 p->numentries=p->numentries+1;
552 p->entries = realloc(p->entries, (p->numentries+1) * sizeof(double) );
553 p->iindices = realloc(p->iindices, (p->numentries+1) * sizeof(int) );
554 p->jindices = realloc(p->jindices, (p->numentries+1) * sizeof(int) );
556 p->entries[p->numentries]=ent;
557 p->iindices[p->numentries]=indexi;
558 p->jindices[p->numentries]=indexj;
567 * If we get here, we have a non-empty structure but not the right block
568 * inside hence create a new structure.
571 p=(struct sparseblock *)calloc(1, sizeof(struct sparseblock));
573 //two entries because this library ignores indices starting in zerox
574 p->constraintnum=matno;
579 p->entries=calloc(p->numentries+1, sizeof(double));
580 p->iindices=calloc(p->numentries+1, sizeof(int));
581 p->jindices=calloc(p->numentries+1, sizeof(int));
583 p->entries[p->numentries]=ent;
584 p->iindices[p->numentries]=indexi;
585 p->jindices[p->numentries]=indexj;
587 p->blocksize=blocksize;
593 int blksz=C->blocks[blkno].blocksize;
594 if (C->blocks[blkno].blockcategory == DIAG){
595 C->blocks[blkno].data.vec[indexi]=ent;
597 C->blocks[blkno].data.mat[ijtok(indexi,indexj,blksz)]=ent;
598 C->blocks[blkno].data.mat[ijtok(indexj,indexi,blksz)]=ent;