DSDP
dsdpadddata.c
Go to the documentation of this file.
1#include "dsdpsdp.h"
2#include "dsdpsys.h"
7#undef __FUNCT__
8#define __FUNCT__ "SDPConeCheckI"
15int SDPConeCheckI(SDPCone sdpcone,int vari){
16 DSDPFunctionBegin;
17 SDPConeValid(sdpcone);
18 if (vari<0 || vari>sdpcone->m) {
19 DSDPSETERR2(1,"Bad Data Matrix: variable: %d (Max: %d)\n",vari,sdpcone->m+1);}
20 DSDPFunctionReturn(0);
21}
22
23#undef __FUNCT__
24#define __FUNCT__ "SDPConeCheckJ"
31int SDPConeCheckJ(SDPCone sdpcone,int blockj){
32 DSDPFunctionBegin;
33 SDPConeValid(sdpcone);
34 if (blockj<0 || blockj>= sdpcone->nblocks) {
35 DSDPSETERR2(2,"Bad Data Matrix: Block: %d (Max: %d)\n",blockj,sdpcone->nblocks-1);}
36 DSDPFunctionReturn(0);
37}
38
39#undef __FUNCT__
40#define __FUNCT__ "SDPConeCheckN"
48int SDPConeCheckN(SDPCone sdpcone,int blockj, int n){
49 int info;
50 DSDPFunctionBegin;
51 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
52 if (sdpcone->blk[blockj].n==0 && n>0){info=SDPConeSetBlockSize(sdpcone,blockj,n);DSDPCHKERR(info);}
53 if (sdpcone->blk[blockj].n != n){
54 DSDPSETERR3(3,"Check Dimension of Data Matrix: Block: %d, %d -- expecting %d\n",
55 blockj,n,sdpcone->blk[blockj].n);
56 }
57 DSDPFunctionReturn(0);
58}
59
60#undef __FUNCT__
61#define __FUNCT__ "SDPConeCheckM"
68int SDPConeCheckM(SDPCone sdpcone,int m){
69 DSDPFunctionBegin;
70 SDPConeValid(sdpcone);
71 if (m!=sdpcone->m){
72 DSDPSETERR1(4,"Check dimension of array. This problem has %d variables\n",sdpcone->m);}
73 DSDPFunctionReturn(0);
74}
75
76#undef __FUNCT__
77#define __FUNCT__ "SDPConeValidStorageFormat"
84int SDPConeValidStorageFormat(SDPCone sdpcone, char format){
85 DSDPFunctionBegin;
86 if (format!='P' && format != 'U'){
87 DSDPSETERR1(4,"Check format of Block: %c is not supported! Use P or U. \n",format);
88 }
89 DSDPFunctionReturn(0);
90}
91
92#undef __FUNCT__
93#define __FUNCT__ "SDPConeCheckStorageFormat"
101int SDPConeCheckStorageFormat(SDPCone sdpcone,int blockj, char format){
102 int info;
103 DSDPFunctionBegin;
104 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
105 info=SDPConeValidStorageFormat(sdpcone,format);DSDPCHKERR(info);
106 if (sdpcone->blk[blockj].format=='N'){
107 sdpcone->blk[blockj].format = format;
108 }
109 if (sdpcone->blk[blockj].format != format){
110 DSDPSETERR3(4,"Check format of Data Matrix: Block: %d, %c -- expecting %c\n",
111 blockj,format,sdpcone->blk[blockj].format);
112 }
113 DSDPFunctionReturn(0);
114}
115
116#undef __FUNCT__
117#define __FUNCT__ "SDPConeRemoveDataMatrix"
127int SDPConeRemoveDataMatrix(SDPCone sdpcone,int blockj, int vari){
128 int info;
129 DSDPFunctionBegin;
130 info=SDPConeCheckI(sdpcone,vari);DSDPCHKERR(info);
131 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
132 info=DSDPBlockRemoveDataMatrix(&sdpcone->blk[blockj].ADATA,vari);DSDPCHKERR(info);
133 DSDPFunctionReturn(0);
134}
135
136#undef __FUNCT__
137#define __FUNCT__ "SDPConeAddDataMatrix"
154int SDPConeAddDataMatrix(SDPCone sdpcone,int blockj, int vari, int n, char format, struct DSDPDataMat_Ops* dsdpdataops, void* data){
155 int info;
156 DSDPFunctionBegin;
157 info=SDPConeCheckI(sdpcone,vari);DSDPCHKERR(info);
158 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
159 info=SDPConeCheckN(sdpcone,blockj,n);DSDPCHKERR(info);
160 info=SDPConeCheckStorageFormat(sdpcone,blockj,format);DSDPCHKERR(info);
161 info=DSDPBlockAddDataMatrix(&sdpcone->blk[blockj].ADATA,vari,dsdpdataops,data);DSDPCHKERR(info);
162 DSDPFunctionReturn(0);
163}
164
165#undef __FUNCT__
166#define __FUNCT__ "SDPConeSetRMatrix"
181int SDPConeSetRMatrix(SDPCone sdpcone,int blockj, int n, char format, struct DSDPDataMat_Ops* dsdpdataops, void* data){
182 int info;
183 int vari=sdpcone->m+1;
184 DSDPFunctionBegin;
185 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
186 info=SDPConeCheckN(sdpcone,blockj,n);DSDPCHKERR(info);
187 info=SDPConeCheckStorageFormat(sdpcone,blockj,format);DSDPCHKERR(info);
188 info=DSDPBlockRemoveDataMatrix(&sdpcone->blk[blockj].ADATA,vari);DSDPCHKERR(info);
189 info=DSDPBlockSetDataMatrix(&sdpcone->blk[blockj].ADATA,vari,dsdpdataops,data);DSDPCHKERR(info);
190 DSDPFunctionReturn(0);
191}
192
193
194#undef __FUNCT__
195#define __FUNCT__ "SDPConeViewDataMatrix"
205int SDPConeViewDataMatrix(SDPCone sdpcone,int blockj, int vari){
206 int info,ii,vari2,nnzmats;
207 DSDPDataMat AA;
208 DSDPFunctionBegin;
209 info=SDPConeCheckI(sdpcone,vari);DSDPCHKERR(info);
210 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
211 info=DSDPBlockCountNonzeroMatrices(&sdpcone->blk[blockj].ADATA,&nnzmats);DSDPCHKERR(info);
212 for (ii=0;ii<nnzmats; ii++){ /* Matrix Entries */
213 info=DSDPBlockGetMatrix(&sdpcone->blk[blockj].ADATA,ii,&vari2,0,&AA);DSDPCHKVARERR(vari,info);
214 if (vari2==vari){ info = DSDPDataMatView(AA);DSDPCHKERR(info);}
215 }
216 DSDPFunctionReturn(0);
217}
218
219#undef __FUNCT__
220#define __FUNCT__ "SDPConeClearVMatrix"
228int SDPConeClearVMatrix(SDPCone sdpcone,int blockj, int n){
229 int info;
230 DSDPFunctionBegin;
231 SDPConeValid(sdpcone);
232 info=DSDPVMatDestroy(&sdpcone->blk[blockj].T);DSDPCHKERR(info);
233 info=DSDPVMatInitialize(&sdpcone->blk[blockj].T);DSDPCHKERR(info);
234 DSDPFunctionReturn(0);
235}
236
237
238#undef __FUNCT__
239#define __FUNCT__ "SDPConeSetXMat"
247int SDPConeSetXMat(SDPCone sdpcone,int blockj, int n){
248 int info;
249 char UPLQ;
250 DSDPVMat T;
251
252 DSDPFunctionBegin;
253 SDPConeValid(sdpcone);
254 info=SDPConeClearVMatrix(sdpcone,blockj,n);DSDPCHKERR(info);
255 DSDPLogInfo(0,10,"Create block X Mat: Block: %d, size: %d.\n",blockj,n);
256 info=SDPConeGetStorageFormat(sdpcone,blockj,&UPLQ); DSDPCHKERR(info);
257 info=DSDPMakeVMat(UPLQ,n,&T);DSDPCHKERR(info);
258 sdpcone->blk[blockj].T=T;
259 DSDPFunctionReturn(0);
260}
261
262#undef __FUNCT__
263#define __FUNCT__ "SDPConeSetXArray"
278int SDPConeSetXArray(SDPCone sdpcone,int blockj, int n, double xx[], int nn){
279 int info;
280 char UPLQ;
281 DSDPVMat T;
282 DSDPFunctionBegin;
283 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
284 info=SDPConeCheckN(sdpcone,blockj,n);DSDPCHKERR(info);
285 info=SDPConeClearVMatrix(sdpcone,blockj,n);DSDPCHKERR(info);
286 DSDPLogInfo(0,10,"Set block X array: Block: %d, size: %d.\n",blockj,n);
287 info=SDPConeGetStorageFormat(sdpcone,blockj,&UPLQ); DSDPCHKERR(info);
288 info=DSDPMakeVMatWithArray(UPLQ,xx,nn,n,&T);DSDPCHKERR(info);
289 sdpcone->blk[blockj].T=T;
290 DSDPFunctionReturn(0);
291}
292
293#undef __FUNCT__
294#define __FUNCT__ "SDPConeGetXArray"
328int SDPConeGetXArray(SDPCone sdpcone,int blockj, double* xx[], int *nn){
329 int info,flag;
330 DSDPFunctionBegin;
331 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
332 info=DSDPVMatExist(sdpcone->blk[blockj].T,&flag);DSDPCHKERR(info);
333 if (flag==0){
334 DSDPSETERR(6,"No X Array available, Call DSDPSetup() or SDPConeSetXArray.\n");}
335 info=DSDPVMatGetArray(sdpcone->blk[blockj].T,xx,nn);DSDPCHKERR(info);
336 DSDPFunctionReturn(0);
337}
338
339#undef __FUNCT__
340#define __FUNCT__ "SDPConeRestoreXArray"
351int SDPConeRestoreXArray(SDPCone sdpcone,int blockj, double* xx[], int *nn){
352 int info,flag;
353 DSDPFunctionBegin;
354 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
355 info=DSDPVMatExist(sdpcone->blk[blockj].T,&flag);DSDPCHKERR(info);
356 if (flag==0){
357 DSDPSETERR(6,"No X Array available, Call DSDPSetup() or SDPConeSetXArray.\n");}
358 info=DSDPVMatRestoreArray(sdpcone->blk[blockj].T,xx,nn);DSDPCHKERR(info);
359 DSDPFunctionReturn(0);
360}
361
362#undef __FUNCT__
363#define __FUNCT__ "SDPConeMatrixView"
372int SDPConeMatrixView(SDPCone sdpcone, int blockj){
373 int info;
374 DSDPFunctionBegin;
375 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
376 info=DSDPVMatView(sdpcone->blk[blockj].T);DSDPCHKERR(info);
377 DSDPFunctionReturn(0);
378}
379
380#undef __FUNCT__
381#define __FUNCT__ "SDPConeUseFullSymmetricFormat"
414int SDPConeUseFullSymmetricFormat(SDPCone sdpcone, int blockj){
415 int info;
416 DSDPFunctionBegin;
417 info=SDPConeSetStorageFormat(sdpcone,blockj,'U');DSDPCHKERR(info);
418 DSDPFunctionReturn(0);
419}
420
421#undef __FUNCT__
422#define __FUNCT__ "SDPConeUsePackedFormat"
452int SDPConeUsePackedFormat(SDPCone sdpcone, int blockj){
453 int info;
454 DSDPFunctionBegin;
455 info=SDPConeSetStorageFormat(sdpcone,blockj,'P');DSDPCHKERR(info);
456 DSDPFunctionReturn(0);
457}
458
459#undef __FUNCT__
460#define __FUNCT__ "SDPConeSetStorageFormat"
479int SDPConeSetStorageFormat(SDPCone sdpcone, int blockj, char format){
480 int info;
481 DSDPFunctionBegin;
482 info=SDPConeValidStorageFormat(sdpcone,format);DSDPCHKERR(info);
483 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
484 sdpcone->blk[blockj].format=format;
485 DSDPFunctionReturn(0);
486}
487
488#undef __FUNCT__
489#define __FUNCT__ "SDPConeGetStorageFormat"
505int SDPConeGetStorageFormat(SDPCone sdpcone, int blockj, char *format){
506 int info;
507 DSDPFunctionBegin;
508 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
509 *format=sdpcone->blk[blockj].format;
510 if (*format=='N') *format='P';
511 DSDPFunctionReturn(0);
512}
513
514#undef __FUNCT__
515#define __FUNCT__ "SDPConeScaleBarrier"
516int SDPConeScaleBarrier(SDPCone sdpcone,int blockj, double ggamma){
517 int info;
518 DSDPFunctionBegin;
519 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
520 sdpcone->blk[blockj].gammamu=ggamma;
521 DSDPFunctionReturn(0);
522}
523
524#undef __FUNCT__
525#define __FUNCT__ "SDPConeSetBlockSize"
535int SDPConeSetBlockSize(SDPCone sdpcone, int blockj, int n){
536 int info,n0;
537 DSDPFunctionBegin;
538 DSDPLogInfo(0,10,"Set block size: Block: %d, size: %d.\n",blockj,n);
539 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
540 n0=sdpcone->blk[blockj].n;
541 if (n0==n){DSDPFunctionReturn(0);}
542 if (n0!=0 &&n0!=n){
543 DSDPSETERR2(5,"Block %d Size previously set to %d \n",blockj,n0); }
544 sdpcone->blk[blockj].n=n;
545 sdpcone->nn+=n-n0;
546 DSDPFunctionReturn(0);
547}
548
549#undef __FUNCT__
550#define __FUNCT__ "SDPConeGetBlockSize"
560int SDPConeGetBlockSize(SDPCone sdpcone, int blockj, int *n){
561 int info;
562 DSDPFunctionBegin;
563 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
564 *n=sdpcone->blk[blockj].n;
565 DSDPFunctionReturn(0);
566}
567
576#undef __FUNCT__
577#define __FUNCT__ "SDPConeGetNumberOfBlocks"
578int SDPConeGetNumberOfBlocks(SDPCone sdpcone, int *nblocks){
579 DSDPFunctionBegin;
580 SDPConeValid(sdpcone);
581 *nblocks=sdpcone->nblocks;
582 DSDPFunctionReturn(0);
583}
584
585#undef __FUNCT__
586#define __FUNCT__ "SDPConeSetSparsity"
596int SDPConeSetSparsity(SDPCone sdpcone, int blockj, int nnz){
597 int info;
598 DSDPFunctionBegin;
599 DSDPLogInfo(0,10,"Set block nonzeros: Block: %d, Nonzero Matrices: %d.\n",blockj,nnz);
600 info=SDPConeCheckJ(sdpcone,blockj);DSDPCHKERR(info);
601 if (nnz>sdpcone->m) nnz=sdpcone->m;
602 info=DSDPBlockDataAllocate(&sdpcone->blk[blockj].ADATA,nnz+2); DSDPCHKERR(info);
603 DSDPFunctionReturn(0);
604}
605
606
607#undef __FUNCT__
608#define __FUNCT__ "SDPConeView"
617int SDPConeView(SDPCone sdpcone){
618 int blockj,info;
619 DSDPFunctionBegin;
620 for (blockj=0; blockj<sdpcone->nblocks; blockj++){
621 printf("Block: %d, Dimension: %d\n",blockj,sdpcone->blk[blockj].n);
622 info=DSDPBlockView(&sdpcone->blk[blockj].ADATA);DSDPCHKERR(info);
623 }
624 DSDPFunctionReturn(0);
625}
626
627#undef __FUNCT__
628#define __FUNCT__ "SDPConeView2"
637int SDPConeView2(SDPCone sdpcone){
638 int blockj,info;
639 DSDPFunctionBegin;
640 for (blockj=0; blockj<sdpcone->nblocks; blockj++){
641 printf("Block: %d, Dimension: %d\n",blockj,sdpcone->blk[blockj].n);
642 info=DSDPBlockView2(&sdpcone->blk[blockj].ADATA);DSDPCHKERR(info);
643 }
644 DSDPFunctionReturn(0);
645}
646
647#undef __FUNCT__
648#define __FUNCT__ "SDPConeView3"
657int SDPConeView3(SDPCone sdpcone){
658 int blockj,id,n,info,nnzmats;
659 DSDPFunctionBegin;
660 for (blockj=0; blockj<sdpcone->nblocks; blockj++){
661 n=sdpcone->blk[blockj].n;
662 printf("Block: %d \n",blockj);
663 printf(" Dimension: %d\n",n);
664 info=DSDPDSMatGetType(sdpcone->blk[blockj].DS,&id);
665 if (id==1){
666 printf(" DS Matrix Type: Dense, Using LAPACK\n");
667 } else {
668 printf(" DS Matrix Type: %d\n",id);
669 }
670 info=DSDPDualMatGetType(sdpcone->blk[blockj].S,&id);
671 if (id==1){
672 printf(" Dual Matrix Type: Dense, Using LAPACK\n");
673 } else {
674 printf(" Dual Matrix Type: %d\n",id);
675 }
676 info=DSDPBlockCountNonzeroMatrices(&sdpcone->blk[blockj].ADATA,&nnzmats);DSDPCHKERR(info);
677 printf(" Number of Data Matrices: %d of %d\n",nnzmats-1,sdpcone->m+1);
678 printf(" Number of Data Nonzeros: %d\n",sdpcone->blk[blockj].nnz);
679 }
680 DSDPFunctionReturn(0);
681}
682
683
684#undef __FUNCT__
685#define __FUNCT__ "SDPConeCheckData"
693 int i,ii,blockj,nnzmats,info;
694 double scl=0;
695 DSDPDataMat AA;
696 DSDPIndex IS;
697 DSDPVMat T;
698 DSDPDSMat DS;
699 DSDPDualMat S1,S2;
700 SDPConeVec W,W2;
701 DSDPFunctionBegin;
702 for (blockj=0; blockj<sdpcone->nblocks; blockj++){
703 T=sdpcone->blk[blockj].T;DS=sdpcone->blk[blockj].DS;
704 W=sdpcone->blk[blockj].W;W2=sdpcone->blk[blockj].W2;
705 S1=sdpcone->blk[blockj].S;S2=sdpcone->blk[blockj].SS;
706 IS=sdpcone->blk[blockj].IS;
707 printf("Block: %d\n",blockj);
708 info=DSDPVMatCheck(T,W,W2);DSDPCHKERR(info);
709 info=DSDPDSMatCheck(DS,W,W2,T);DSDPCHKERR(info);
710 info=DSDPDualMatCheck(S1,W,W2,IS,T);DSDPCHKERR(info);
711 info=DSDPDualMatCheck(S2,W,W2,IS,T);DSDPCHKERR(info);
712
713 info=DSDPBlockCountNonzeroMatrices(&sdpcone->blk[blockj].ADATA,&nnzmats);DSDPCHKERR(info);
714 for (ii=0;ii<nnzmats;ii++){
715 info=DSDPBlockGetMatrix(&sdpcone->blk[blockj].ADATA,ii,&i,&scl,&AA);DSDPCHKERR(info);
716 if (i==0) continue;
717 printf(" Variable: %d, \n",i);
718 info=DSDPDataMatCheck(AA,W,IS,T);DSDPCHKERR(info);
719 }
720 }
721 DSDPFunctionReturn(0);
722}
struct SDPCone_C * SDPCone
The SDPCone object points to blocks of data that specify semidefinite matrix inequalities.
Definition dsdp5.h:26
int SDPConeView(SDPCone sdpcone)
Print the SDPCone to the screen;.
int SDPConeClearVMatrix(SDPCone sdpcone, int blockj, int n)
Free V matrix.
int SDPConeSetRMatrix(SDPCone sdpcone, int blockj, int n, char format, struct DSDPDataMat_Ops *dsdpdataops, void *data)
Add identity to dual matrix.
int SDPConeSetXMat(SDPCone sdpcone, int blockj, int n)
Create X matrix.
int SDPConeCheckJ(SDPCone sdpcone, int blockj)
Check validity of parameter.
Definition dsdpadddata.c:31
int SDPConeCheckI(SDPCone sdpcone, int vari)
Check validity of parameter.
Definition dsdpadddata.c:15
int SDPConeCheckN(SDPCone sdpcone, int blockj, int n)
Check validity of parameter.
Definition dsdpadddata.c:48
int SDPConeCheckM(SDPCone sdpcone, int m)
Check validity of parameter.
Definition dsdpadddata.c:68
int SDPConeValidStorageFormat(SDPCone sdpcone, char format)
Check validity of parameter.
Definition dsdpadddata.c:84
int SDPConeCheckStorageFormat(SDPCone sdpcone, int blockj, char format)
Check validity of parameters.
int DSDPDataMatCheck(DSDPDataMat AA, SDPConeVec W, DSDPIndex IS, DSDPVMat XX)
Check correctness of operations on the data.
Definition dsdpblock.c:498
int DSDPBlockView2(DSDPBlockData *ADATA)
Print the data.
Definition dsdpblock.c:474
int DSDPBlockDataAllocate(DSDPBlockData *ADATA, int nnz)
Allocate some structures.
Definition dsdpblock.c:221
int DSDPBlockGetMatrix(DSDPBlockData *ADATA, int id, int *vari, double *scl, DSDPDataMat *A)
Get a data matrix from a block of data.
Definition dsdpblock.c:307
int DSDPBlockAddDataMatrix(DSDPBlockData *ADATA, int vari, struct DSDPDataMat_Ops *dsdpdataops, void *data)
Add data matrix into SDP block.
Definition dsdpblock.c:381
int DSDPBlockSetDataMatrix(DSDPBlockData *ADATA, int vari, struct DSDPDataMat_Ops *dsdpdataops, void *data)
Set data matrix into SDP block.
Definition dsdpblock.c:406
int DSDPBlockRemoveDataMatrix(DSDPBlockData *ADATA, int vari)
Remove a data matrix.
Definition dsdpblock.c:351
int DSDPBlockCountNonzeroMatrices(DSDPBlockData *ADATA, int *nzmats)
Count how many data matrices are in a block of data.
Definition dsdpblock.c:272
int DSDPBlockView(DSDPBlockData *ADATA)
Print the structure of the block.
Definition dsdpblock.c:454
int DSDPDataMatView(DSDPDataMat A)
Print matrix.
struct DSDPDataMat_C DSDPDataMat
Represents a single symmetric data matrix for one block in this semidefinite cone.
Definition dsdpdatamat.h:25
struct DSDPDSMat_C DSDPDSMat
A symmetric Delta S matrix for one block in the semidefinite cone.
Definition dsdpdsmat.h:33
struct DSDPDualMat_C DSDPDualMat
Represents an S matrix for one block in the semidefinite cone.
Definition dsdpdualmat.h:27
Internal SDPCone data structures and routines.
int DSDPMakeVMatWithArray(char, double[], int, int, DSDPVMat *)
Allocate V matrix using the given array.
Definition sdpsss.c:381
int DSDPMakeVMat(char, int, DSDPVMat *)
Allocate V matrix.
Definition sdpsss.c:351
Error handling, printing, and profiling.
int DSDPVMatInitialize(DSDPVMat *B)
Set pointers to null.
Definition dsdpxmat.c:424
int DSDPVMatGetArray(DSDPVMat X, double **v, int *nn)
Get the array that stores the matrix.
Definition dsdpxmat.c:211
int DSDPVMatRestoreArray(DSDPVMat X, double **v, int *nn)
Restore the array that stores the matrix.
Definition dsdpxmat.c:233
int DSDPVMatExist(DSDPVMat X, int *flag)
Answer whether the array has been allocated or not.
Definition dsdpxmat.c:440
int DSDPVMatCheck(DSDPVMat X, SDPConeVec W1, SDPConeVec W2)
Test correctness of operations.
Definition dsdpxmat.c:327
int DSDPVMatDestroy(DSDPVMat *X)
Deallocate matrix.
Definition dsdpxmat.c:86
int DSDPVMatView(DSDPVMat X)
Print matrix.
Definition dsdpxmat.c:107
struct DSDPVMat_C DSDPVMat
Represents a dense symmetric matrix for one block in the semidefinite cone.
Definition dsdpxmat.h:26
int SDPConeViewDataMatrix(SDPCone sdpcone, int blockj, int vari)
Print a data matrix to the screen.
int SDPConeGetXArray(SDPCone sdpcone, int blockj, double *xx[], int *nn)
After applying the solver, set a pointer to the array in the object with the solution X.
int SDPConeAddDataMatrix(SDPCone sdpcone, int blockj, int vari, int n, char format, struct DSDPDataMat_Ops *dsdpdataops, void *data)
Add a data matrix .
int SDPConeUseFullSymmetricFormat(SDPCone sdpcone, int blockj)
Use full symmetric format for the dense array.
int SDPConeRemoveDataMatrix(SDPCone sdpcone, int blockj, int vari)
Remove the data matrix from the cone.
int SDPConeCheckData(SDPCone sdpcone)
Check the matrix operations on a data matrix;.
int SDPConeUsePackedFormat(SDPCone sdpcone, int blockj)
Use packed symmetric format for the dense array.
int SDPConeView3(SDPCone sdpcone)
Print the SDP cone to the screen in a third way.
int SDPConeView2(SDPCone sdpcone)
Print the SDP cone to the screen in a second way.
int SDPConeGetStorageFormat(SDPCone sdpcone, int blockj, char *format)
Get the storage format for the block.
int SDPConeSetStorageFormat(SDPCone sdpcone, int blockj, char format)
Set the dense storage format of a block in the semidefinite cone.
int SDPConeSetBlockSize(SDPCone sdpcone, int blockj, int n)
Set the dimension of one block in the semidefinite cone.
int SDPConeGetNumberOfBlocks(SDPCone sdpcone, int *nblocks)
Get the number of blocks in the semidefinite cone.
int SDPConeSetSparsity(SDPCone sdpcone, int blockj, int nnz)
Set the number of nonzero matrices in a block of the semidefinite cone.
int SDPConeSetXArray(SDPCone sdpcone, int blockj, int n, double xx[], int nn)
Provide an array for the SDPCone object can use to store dense matrices.
int SDPConeMatrixView(SDPCone sdpcone, int blockj)
Print the dense array to the screen.
int SDPConeGetBlockSize(SDPCone sdpcone, int blockj, int *n)
Get the dimension of one block in the semidefinite cone.
int SDPConeRestoreXArray(SDPCone sdpcone, int blockj, double *xx[], int *nn)
Restore the dense array and set these pointers to null.
struct SDPConeVec_C SDPConeVec
SDPConeVec is a vector with the dimension of the block in the SDP cone.
Definition sdpconevec.h:26
Table of function pointers that operate on the data matrix.