Limbo 3.5.4
Loading...
Searching...
No Matches
CsdpEasySdpApi.h
Go to the documentation of this file.
1
10/*
11 * This is an easy to call version of the sdp routine. It takes as
12 * input a problem (n,k,C,a,constraints,constant_offset), and an
13 * initial solution (X,y,Z), allocates working storage, and calls sdp()
14 * to solve the problem. The solution is returned in X,y,Z,pobj,dobj, and
15 * the return code from sdp is returned as the return value from easy_sdp.
16 *
17 */
18
19#ifndef LIMBO_SOLVERS_CSDPEASYSDPAPI_H
20#define LIMBO_SOLVERS_CSDPEASYSDPAPI_H
21
22// must define NOSHORTS to forbid the usage of unsigned shorts in Csdp
23#define NOSHORTS
24
25// necessary since we are linking C library
26#ifdef __cplusplus
27extern "C" {
28#endif
29#include <stdlib.h>
30#include <stdio.h>
31#include <math.h>
32#include <limbo/thirdparty/Csdp/include/declarations.h>
33#ifdef __cplusplus
34}
35#endif
36
38namespace limbo
39{
41namespace solvers
42{
43
58template <typename T>
60 // this part is same as original version
61 int n,
62 int k,
63 struct blockmatrix C,
64 double *a,
65 struct constraintmatrix *constraints,
66 double constant_offset,
67 struct blockmatrix *pX,
68 double **py,
69 struct blockmatrix *pZ,
70 double *ppobj,
71 double *pdobj,
72 // newly added parameters
73 struct paramstruc const& params,
74 int const& printlevel
75 )
76{
77 int ret;
78 struct constraintmatrix fill;
79 //struct paramstruc params; // commence out because we pass it as a parameter
80 struct blockmatrix work1;
81 struct blockmatrix work2;
82 struct blockmatrix work3;
83 struct blockmatrix bestx;
84 struct blockmatrix bestz;
85 struct blockmatrix Zi;
86 struct blockmatrix dZ;
87 struct blockmatrix dX;
88 struct blockmatrix cholxinv;
89 struct blockmatrix cholzinv;
90 double *workvec1;
91 double *workvec2;
92 double *workvec3;
93 double *workvec4;
94 double *workvec5;
95 double *workvec6;
96 double *workvec7;
97 double *workvec8;
98 double *diagO;
99 double *Fp;
100 double *O;
101 double *dy;
102 double *dy1;
103 double *rhs;
104 double *besty;
105 //int printlevel; // commence out because we pass it as a parameter
106 int ldam;
107 struct sparseblock **byblocks;
108 struct sparseblock *ptr;
109 struct sparseblock *oldptr;
110 int i;
111 int j;
112 //int blk; // commence out because compiler reports set but unused variable
113 struct sparseblock *p;
114 struct sparseblock *q;
115 struct sparseblock *prev=NULL;
116 double gap;
117 int nnz;
118
119 /*
120 * Initialize the parameters.
121 */
122
123 //initparams(&params,&printlevel); // commence out because we pass it as a parameter
124
125 /*
126 * Allocate working storage
127 */
128
129 alloc_mat(C,&work1);
130 alloc_mat(C,&work2);
131 alloc_mat(C,&work3);
132 alloc_mat_packed(C,&bestx);
133 alloc_mat_packed(C,&bestz);
134 alloc_mat_packed(C,&cholxinv);
135 alloc_mat_packed(C,&cholzinv);
136
137 besty=(double *)malloc(sizeof(double)*(k+1));
138 if (besty == NULL)
139 {
140 printf("Storage Allocation Failed!\n");
141 exit(10);
142 };
143
144 if (n > k)
145 {
146 workvec1=(double *)malloc(sizeof(double)*(n+1));
147 }
148 else
149 {
150 workvec1=(double *)malloc(sizeof(double)*(k+1));
151 };
152 if (workvec1 == NULL)
153 {
154 printf("Storage Allocation Failed!\n");
155 exit(10);
156 };
157
158 if (n > k)
159 {
160 workvec2=(double *)malloc(sizeof(double)*(n+1));
161 }
162 else
163 {
164 workvec2=(double *)malloc(sizeof(double)*(k+1));
165 };
166 if (workvec2 == NULL)
167 {
168 printf("Storage Allocation Failed!\n");
169 exit(10);
170 };
171
172 if (n > k)
173 {
174 workvec3=(double *)malloc(sizeof(double)*(n+1));
175 }
176 else
177 {
178 workvec3=(double *)malloc(sizeof(double)*(k+1));
179 };
180 if (workvec3 == NULL)
181 {
182 printf("Storage Allocation Failed!\n");
183 exit(10);
184 };
185
186 if (n > k)
187 {
188 workvec4=(double *)malloc(sizeof(double)*(n+1));
189 }
190 else
191 {
192 workvec4=(double *)malloc(sizeof(double)*(k+1));
193 };
194 if (workvec4 == NULL)
195 {
196 printf("Storage Allocation Failed!\n");
197 exit(10);
198 };
199
200 if (n > k)
201 {
202 workvec5=(double *)malloc(sizeof(double)*(n+1));
203 }
204 else
205 {
206 workvec5=(double *)malloc(sizeof(double)*(k+1));
207 };
208 if (workvec5 == NULL)
209 {
210 printf("Storage Allocation Failed!\n");
211 exit(10);
212 };
213
214 if (n > k)
215 {
216 workvec6=(double *)malloc(sizeof(double)*(n+1));
217 }
218 else
219 {
220 workvec6=(double *)malloc(sizeof(double)*(k+1));
221 };
222 if (workvec6 == NULL)
223 {
224 printf("Storage Allocation Failed!\n");
225 exit(10);
226 };
227
228 if (n > k)
229 {
230 workvec7=(double *)malloc(sizeof(double)*(n+1));
231 }
232 else
233 {
234 workvec7=(double *)malloc(sizeof(double)*(k+1));
235 };
236 if (workvec7 == NULL)
237 {
238 printf("Storage Allocation Failed!\n");
239 exit(10);
240 };
241
242 if (n > k)
243 {
244 workvec8=(double *)malloc(sizeof(double)*(n+1));
245 }
246 else
247 {
248 workvec8=(double *)malloc(sizeof(double)*(k+1));
249 };
250 if (workvec8 == NULL)
251 {
252 printf("Storage Allocation Failed!\n");
253 exit(10);
254 };
255
256
257 if (n > k)
258 {
259 diagO=(double *)malloc(sizeof(double)*(n+1));
260 }
261 else
262 {
263 diagO=(double *)malloc(sizeof(double)*(k+1));
264 };
265 if (diagO == NULL)
266 {
267 printf("Storage Allocation Failed!\n");
268 exit(10);
269 };
270
271
272
273 rhs=(double*)malloc(sizeof(double)*(k+1));
274 if (rhs == NULL)
275 {
276 printf("Storage Allocation Failed!\n");
277 exit(10);
278 };
279
280 dy=(double*)malloc(sizeof(double)*(k+1));
281 if (dy == NULL)
282 {
283 printf("Storage Allocation Failed!\n");
284 exit(10);
285 };
286
287 dy1=(double*)malloc(sizeof(double)*(k+1));
288 if (dy1 == NULL)
289 {
290 printf("Storage Allocation Failed!\n");
291 exit(10);
292 };
293
294 Fp=(double*)malloc(sizeof(double)*(k+1));
295 if (Fp == NULL)
296 {
297 printf("Storage Allocation Failed!\n");
298 exit(10);
299 };
300
301 /*
302 * Work out the leading dimension for the array. Note that we may not
303 * want to use k itself, for cache issues.
304 */
305 if ((k % 2) == 0)
306 ldam=k+1;
307 else
308 ldam=k;
309
310 O=(double*)malloc(sizeof(double)*ldam*ldam);
311 if (O == NULL)
312 {
313 printf("Storage Allocation Failed!\n");
314 exit(10);
315 };
316
317 alloc_mat(C,&Zi);
318 alloc_mat(C,&dZ);
319 alloc_mat(C,&dX);
320
321 /*
322 * Fill in lots of details in the constraints data structure that haven't
323 * necessarily been done before now.
324 */
325
326 /*
327 * Set up the cross links used by op_o
328 * While we're at it, determine which blocks are sparse and dense.
329 */
330
331 /*
332 * Next, setup issparse and NULL out all nextbyblock pointers.
333 */
334
335 for (i=1; i<=k; i++)
336 {
337 p=constraints[i].blocks;
338 while (p != NULL)
339 {
340 /*
341 * First, set issparse.
342 */
343 if (((p->numentries) > 0.25*(p->blocksize)) && ((p->numentries) > 15))
344 {
345 p->issparse=0;
346 }
347 else
348 {
349 p->issparse=1;
350 };
351
352 if (C.blocks[p->blocknum].blockcategory == DIAG)
353 p->issparse=1;
354
355 /*
356 * Setup the cross links.
357 */
358
359 p->nextbyblock=NULL;
360 p=p->next;
361 };
362 };
363
364 /*
365 * Now, cross link.
366 */
367 for (i=1; i<=k; i++)
368 {
369 p=constraints[i].blocks;
370 while (p != NULL)
371 {
372 if (p->nextbyblock == NULL)
373 {
374 //blk=p->blocknum; // commence out because compiler reports set but unused variable
375
376 /*
377 * link in the remaining blocks.
378 */
379 for (j=i+1; j<=k; j++)
380 {
381 q=constraints[j].blocks;
382
383 while (q != NULL)
384 {
385 if (q->blocknum == p->blocknum)
386 {
387 if (p->nextbyblock == NULL)
388 {
389 p->nextbyblock=q;
390 q->nextbyblock=NULL;
391 prev=q;
392 }
393 else
394 {
395 prev->nextbyblock=q;
396 q->nextbyblock=NULL;
397 prev=q;
398 };
399 break;
400 };
401 q=q->next;
402 };
403 };
404 };
405 p=p->next;
406 };
407 };
408
409 /*
410 * If necessary, print out information on sparsity of blocks.
411 */
412
413 if (printlevel >= 4)
414 {
415 for (i=1; i<=k; i++)
416 {
417 p=constraints[i].blocks;
418 while (p != NULL)
419 {
420 printf("%d,%d,%d,%d \n",i,p->blocknum,p->issparse,p->numentries);
421 p=p->next;
422 };
423 };
424 };
425
426 /*
427 * Allocate space for byblocks pointers.
428 */
429
430 byblocks=(struct sparseblock **)malloc((C.nblocks+1)*sizeof(struct sparseblock *));
431 if (byblocks == NULL)
432 {
433 printf("Storage Allocation Failed!\n");
434 exit(10);
435 };
436
437 for (i=1; i<=C.nblocks; i++)
438 byblocks[i]=NULL;
439
440 /*
441 * Fill in byblocks pointers.
442 */
443 for (i=1; i<=k; i++)
444 {
445 ptr=constraints[i].blocks;
446 while (ptr != NULL)
447 {
448 if (byblocks[ptr->blocknum]==NULL)
449 {
450 byblocks[ptr->blocknum]=ptr;
451 };
452 ptr=ptr->next;
453 };
454 };
455
456 /*
457 * Compute "fill". This data structure tells us which elements in the
458 * block diagonal matrix have corresponding elements in one of the
459 * constraints, and which constraint this element first appears in.
460 *
461 */
462
463 makefill(k,C,constraints,&fill,work1,printlevel);
464
465 /*
466 * Compute the nonzero structure of O.
467 */
468
469 nnz=structnnz(n,k,C,constraints);
470
471 if (printlevel >= 3)
472 printf("Structural density of O %d, %e \n",nnz,nnz*1.0/(k*k*1.0));
473
474 /*
475 * Sort entries in diagonal blocks of constraints.
476 */
477
478 sort_entries(k,C,constraints);
479
480 /*
481 * Now, call sdp().
482 */
483
484 ret=sdp(n,k,C,a,constant_offset,constraints,byblocks,fill,*pX,*py,*pZ,
485 cholxinv,cholzinv,ppobj,pdobj,work1,work2,work3,workvec1,
486 workvec2,workvec3,workvec4,workvec5,workvec6,workvec7,workvec8,
487 diagO,bestx,besty,bestz,Zi,O,rhs,dZ,dX,dy,dy1,Fp,
488 printlevel,params);
489
490 if (printlevel >= 1)
491 {
492 if (ret==0)
493 printf("Success: SDP solved\n");
494 if (ret==1)
495 printf("Success: SDP is primal infeasible\n");
496 if (ret==2)
497 printf("Success: SDP is dual infeasible\n");
498 if (ret==3)
499 printf("Partial Success: SDP solved with reduced accuracy\n");
500 if (ret >= 4)
501 printf("Failure: return code is %d \n",ret);
502
503 if (ret==1)
504 {
505 op_at(k,*py,constraints,work1);
506 addscaledmat(work1,-1.0,*pZ,work1);
507 printf("Certificate of primal infeasibility: a'*y=%.5e, ||A'(y)-Z||=%.5e\n",-1.0,Fnorm(work1));
508 };
509
510 if (ret==2)
511 {
512 op_a(k,constraints,*pX,workvec1);
513 printf("Certificate of dual infeasibility: tr(CX)=%.5e, ||A(X)||=%.5e\n",trace_prod(C,*pX),norm2(k,workvec1+1));
514 };
515
516 if ((ret==0) || (ret>=3))
517 {
518 if (printlevel >= 3)
519 {
520 printf("XZ Gap: %.7e \n",trace_prod(*pZ,*pX));
521 gap=*pdobj-*ppobj;
522 printf("Real Gap: %.7e \n",gap);
523 };
524
525 if (printlevel >= 1)
526 {
527 gap=*pdobj-*ppobj;
528
529 printf("Primal objective value: %.7e \n",*ppobj);
530 printf("Dual objective value: %.7e \n",*pdobj);
531 printf("Relative primal infeasibility: %.2e \n",
532 pinfeas(k,constraints,*pX,a,workvec1));
533 printf("Relative dual infeasibility: %.2e \n",
534 dinfeas(k,C,constraints,*py,*pZ,work1));
535 printf("Real Relative Gap: %.2e \n",gap/(1+fabs(*pdobj)+fabs(*ppobj)));
536 printf("XZ Relative Gap: %.2e \n",trace_prod(*pZ,*pX)/(1+fabs(*pdobj)+fabs(*ppobj)));
537
538 printf("DIMACS error measures: %.2e %.2e %.2e %.2e %.2e %.2e\n",
539 pinfeas(k,constraints,*pX,a,workvec1)*(1+norm2(k,a+1))/
540 (1+norminf(k,a+1)),
541 0.0,
542 dimacserr3(k,C,constraints,*py,*pZ,work1),
543 0.0,
544 gap/(1+fabs(*pdobj)+fabs(*ppobj)),
545 trace_prod(*pZ,*pX)/(1+fabs(*pdobj)+fabs(*ppobj)));
546 };
547 };
548
549 };
550
551 /*
552 * Now, free up all of the storage.
553 */
554
555 free_mat(work1);
556 free_mat(work2);
557 free_mat(work3);
558 free_mat_packed(bestx);
559 free_mat_packed(bestz);
560 free_mat_packed(cholxinv);
561 free_mat_packed(cholzinv);
562
563 free_mat(Zi);
564 free_mat(dZ);
565 free_mat(dX);
566
567 free(besty);
568 free(workvec1);
569 free(workvec2);
570 free(workvec3);
571 free(workvec4);
572 free(workvec5);
573 free(workvec6);
574 free(workvec7);
575 free(workvec8);
576 free(rhs);
577 free(dy);
578 free(dy1);
579 free(Fp);
580 free(O);
581 free(diagO);
582 free(byblocks);
583
584 /*
585 * Free up the fill data structure.
586 */
587
588 ptr=fill.blocks;
589 while (ptr != NULL)
590 {
591 free(ptr->entries);
592 free(ptr->iindices);
593 free(ptr->jindices);
594 oldptr=ptr;
595 ptr=ptr->next;
596 free(oldptr);
597 };
598
599
600 /*
601 * Finally, free the constraints array.
602 */
603
604 return(ret);
605
606}
607
608}} // namespace limbo // namespace solvers
609
610#endif
mathematical utilities such as abs
namespace for Limbo.Solvers
int easy_sdp_ext(int n, int k, struct blockmatrix C, double *a, struct constraintmatrix *constraints, double constant_offset, struct blockmatrix *pX, double **py, struct blockmatrix *pZ, double *ppobj, double *pdobj, struct paramstruc const &params, int const &printlevel)
API to call Csdp solver.
namespace for Limbo