File: declarations.h

package info (click to toggle)
coinor-csdp 6.2.0-6
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 2,048 kB
  • sloc: ansic: 7,195; makefile: 99
file content (293 lines) | stat: -rw-r--r-- 10,022 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
/*
 * CSDPDECLARATIONS is used to prevent redefinitions if this file is included
 * twice.
 */

#ifndef CSDPDECLARATIONS
#define CSDPDECLARATIONS 
#ifdef __cplusplus
extern "C" {
#endif
/*
  Other important includes that we need.
 */

#include <csdp/index.h>
#include <csdp/blockmat.h>
#include <csdp/parameters.h>

/*
  Our own routines.
  */

void triu(struct blockmatrix A);

void store_packed(struct blockmatrix A, struct blockmatrix B);

void store_unpacked(struct blockmatrix A, struct blockmatrix B);

void alloc_mat_packed(struct blockmatrix A, struct blockmatrix *pB);

void free_mat_packed(struct blockmatrix A);

int structnnz(int n, int k, struct blockmatrix C, const struct constraintmatrix *constraints);

int actnnz(int n, int lda, const double *A);

int bandwidth(int n, int lda, const double *A);

void qreig(int n, double *maindiag, double *offdiag);

void sort_entries(int k, struct blockmatrix C, const struct constraintmatrix *constraints);

double norm2(int n, double *x);

double norm1(int n, double *x);

double norminf(int n, double *x);

double Fnorm(struct blockmatrix A);

double Knorm(struct blockmatrix A);

double mat1norm(struct blockmatrix A);

double matinfnorm(struct blockmatrix A);

double calc_pobj(struct blockmatrix C, struct blockmatrix X, 
		 double constant_offset);

double calc_dobj(int k, double *a, double *y, double constant_offset);

double trace_prod(struct blockmatrix A, struct blockmatrix B);

double linesearch(int n, struct blockmatrix dX,
		  struct blockmatrix work1, struct blockmatrix work2, 
		  struct blockmatrix work3, struct blockmatrix cholinv, 
		  double *q, double *z, double *workvec,
		  double stepfrac,double start, int printlevel);

double pinfeas(int k, const struct constraintmatrix *constraints,
	       struct blockmatrix X, double *a, double *workvec);

double dinfeas(int k, struct blockmatrix C, 
	       const struct constraintmatrix *constraints, const double *y, 
	       struct blockmatrix Z, struct blockmatrix work1);

double dimacserr3(int k, struct blockmatrix C, 
	       const struct constraintmatrix *constraints, const double *y, 
	       struct blockmatrix Z, struct blockmatrix work1);

void op_a(int k, const struct constraintmatrix *constraints,
	  struct blockmatrix X, double *result);

void op_at(int k, const double *y, const struct constraintmatrix *constraints,
	   struct blockmatrix result);

void makefill(int k, struct blockmatrix C, 
	      const struct constraintmatrix *constraints, 
	      struct constraintmatrix *pfill, struct blockmatrix work1, 
	      int printlevel);

void op_o(int k, const struct constraintmatrix *constraints,
	  const struct sparseblock **byblocks, struct blockmatrix Zi, 
          struct blockmatrix X, double *O, struct blockmatrix work1, 
          struct blockmatrix work2);

void addscaledmat(struct blockmatrix A, double scale, struct blockmatrix B,
		  struct blockmatrix C);

void zero_mat(struct blockmatrix A);

void add_mat(struct blockmatrix A,struct blockmatrix B);

void sym_mat(struct blockmatrix A);

void make_i(struct blockmatrix A);

void copy_mat(struct blockmatrix A, struct blockmatrix B);

void mat_mult(double scale1, double scale2, struct blockmatrix A,
	      struct blockmatrix B, struct blockmatrix C);

void mat_multspa(double scale1, double scale2, struct blockmatrix A,
		 struct blockmatrix B, struct blockmatrix C, 
		 struct constraintmatrix fill);

void mat_multspb(double scale1, double scale2, struct blockmatrix A,
		 struct blockmatrix B, struct blockmatrix C, 
		 struct constraintmatrix fill);

void mat_multspc(double scale1, double scale2, struct blockmatrix A,
		 struct blockmatrix B, struct blockmatrix C, 
		 struct constraintmatrix fill);

void mat_mult_raw(int n, double scale1, double scale2, double *ap,
		  double *bp, double *cp);

void mat_mult_rawatlas(int n, double scale1, double scale2, double *ap,
		  double *bp, double *cp);

void matvec(struct blockmatrix A, double *x, double *y);

void alloc_mat(struct blockmatrix A, struct blockmatrix *pB);

void free_mat(struct blockmatrix A);

void initparams(struct paramstruc *params, int *pprintlevel);

void initsoln(int n, int k, struct blockmatrix C, double *a, 
	      const struct constraintmatrix *constraints, struct blockmatrix *pX0,
	      double **py0, struct blockmatrix *pZ0);

void trans(struct blockmatrix A);

void chol_inv(struct blockmatrix A, struct blockmatrix B);

int chol(struct blockmatrix A);

int solvesys(int m, int ldam, double *A, double *rhs);

int user_exit(int n, int k, struct blockmatrix C, const double *a, double dobj,
	      double pobj, double constant_offset, 
	      const struct constraintmatrix *constraints, struct blockmatrix X,
	      const double *y, struct blockmatrix Z, struct paramstruc params);

int read_sol(const char *fname, int n, int k, struct blockmatrix C, 
	     struct blockmatrix *pX, double **py, struct blockmatrix *pZ);

int read_prob(const char *fname, int *pn, int *pk, struct blockmatrix *pC,
	      double **pa, struct constraintmatrix **pconstraints,
	      int printlevel);

int write_prob(const char *fname, int n, int k, struct blockmatrix C,
	       const double *a, const struct constraintmatrix *constraints);

int write_sol(const char *fname, int n, int k, struct blockmatrix X,
	      const double *y, struct blockmatrix Z);

void free_prob(int n, int k, struct blockmatrix C, double *a, 
	       struct constraintmatrix *constraints, struct blockmatrix X,
	       double *y, struct blockmatrix Z);

int sdp(int n, int k, struct blockmatrix C, double *a, double constant_offset,
	const struct constraintmatrix *constraints, const struct sparseblock **byblocks,
	struct constraintmatrix fill, struct blockmatrix X, double *y, 
	struct blockmatrix Z, struct blockmatrix cholxinv,
	struct blockmatrix cholzinv, double *pobj, double *dobj, 
	struct blockmatrix work1, struct blockmatrix work2,     
	struct blockmatrix work3,
	double *workvec1, double *workvec2,
	double *workvec3, double *workvec4, double *workvec5,
	double *workvec6, const double *workvec7, double *workvec8, 
	double *diagO, struct blockmatrix bestx, double *besty,
	struct blockmatrix bestz, struct blockmatrix Zi, double *O,
	double *rhs, struct blockmatrix dZ, struct blockmatrix dX, 
	double *dy, double *dy1, double *Fp, int printlevel, 
	struct paramstruc parameters);

int easy_sdp(int n, int k, struct blockmatrix C, double *a, 
	     const struct constraintmatrix *constraints, double constant_offset,
	     const struct blockmatrix *pX, double **py, const struct blockmatrix *pZ,
	     double *ppobj, double *pdobj);

int checkconstraints(int n, int k, struct blockmatrix C,
		     const struct constraintmatrix *constraints, int printlevel);

int checkc(int n, struct blockmatrix C, int printlevel);
  
void tweakgap(int n, int k, double *a, const struct constraintmatrix *constraints,
	      double gap, struct blockmatrix Z, struct blockmatrix dZ, 
	      double *y, double *dy, struct blockmatrix work1, 
	      struct blockmatrix work2, struct blockmatrix work3, 
	      struct blockmatrix work4, double *workvec1, double *workvec2,
	      double *workvec3, double *workvec4, int printlevel);

int bisect_(int *n, double *eps1, double *d, double *e, double *e2,
	    double *lb, double *ub, int *mm, int *m, double *w, int *ind, 
	    int *ierr, double *rv4, double *rv5);

/*
  BLAS and LINPACK stuff.
  */

/*
  First, BLAS routines.
 */

#ifdef CAPSBLAS
#ifdef NOUNDERBLAS
#  define BLAS_FUNC(name, NAME) NAME
#else
#  define BLAS_FUNC(name, NAME) NAME##_
#endif
#else
#ifdef NOUNDERBLAS
#  define BLAS_FUNC(name, NAME) name
#else
#  define BLAS_FUNC(name, NAME) name##_
#endif
#endif

double BLAS_FUNC(dnrm2, DNRM2)(const int *n, const double *x, const int *incx);
double BLAS_FUNC(dasum, DASUM)(const int *n, const double *x, const int *incx);
double BLAS_FUNC(ddot, DDOT)(const int *n, const double *x, const int *incx,
			     const double *y, const int *incy);
int BLAS_FUNC(idamax, IDAMAX)(const int *n, const double *x, const int *incx);
void BLAS_FUNC(dgemm, DGEMM)(const char *transa, const char *transb,
			     const int *m, const int *n, const int *k,
			     const double *alpha, const double *a,
			     const int *lda, const double *b, const int *ldb,
			     const double *beta, double *c, const int *ldc);
void BLAS_FUNC(dgemv, DGEMV)(const char *trans, const int *m, const int *n,
			     const double *alpha, const double *a,
			     const int *lda, const double *x, const int *incx,
			     const double *beta, double *y, const int *incy);
void BLAS_FUNC(dger, DGER)(const int *m, const int *n, const double *alpha,
			   const double *x,const int *incx,
			   const double *y, const int *incy,
			   double *a, const int *lda);
void BLAS_FUNC(dtrsm, DTRSM)(const char *side, const char *uplo,
			     const char *transa, const char *diag,
			     const int *m, const int *n,
			     const double *alpha,
			     const double *a, const int *lda,
			     double *b, const int *ldb);
void BLAS_FUNC(dtrmv, DTRMV)(const char *uplo, const char *trans,
			     const char *diag, const int *n, const double *a,
			     const int *lda, double *x, const int *incx);

/*
  LAPACK next.
 */

#ifdef CAPSLAPACK
#ifdef NOUNDERLAPACK
#  define LAPACK_FUNC(name, NAME) NAME
#else
#  define LAPACK_FUNC(name, NAME) NAME##_
#endif
#else
#ifdef NOUNDERLAPACK
#  define LAPACK_FUNC(name, NAME) name
#else
#  define LAPACK_FUNC(name, NAME) name##_
#endif
#endif

void LAPACK_FUNC(dpotrf, DPOTRF)(const char *uplo, const int *n,
				 double *a, const int *lda, int *info);
void LAPACK_FUNC(dpotrs, DPOTRS)(const char *uplo, const int *n,
				 const int *nrhs, const double *a,
				 const int *lda, double *b, const int *ldb,
				 int *info);
void LAPACK_FUNC(dpotri, DPOTRI)(const char *uplo, const int *n,
				 double *a, const int *lda, int *info);
void LAPACK_FUNC(dtrtri, DTRTRI)(const char *uplo, const char *diag,
				 const int *n, double *a, const int *lda,
				 int *info);

#ifdef __cplusplus
}
#endif
#endif