49 #include "blockmemshell/memory.h"
50 #include "scip/type_retcode.h"
56 #define SDPA_VERSION 738
59 #define BMS_CALL(x) do \
63 SCIPerrorMessage("No memory in function call\n"); \
64 return SCIP_NOMEMORY; \
70 #define SCIP_RealTOINT(x) ((LAPACKINTTYPE) (x + 0.5))
81 char* JOBZ,
char* RANGE,
char* UPLO,
83 SCIP_Real* VL, SCIP_Real* VU,
85 SCIP_Real* ABSTOL,
LAPACKINTTYPE* M, SCIP_Real* W, SCIP_Real* Z,
115 SCIP_Bool geteigenvectors,
119 SCIP_Real* eigenvalue,
120 SCIP_Real* eigenvector
124 #if ( SDPA_VERSION == 740 )
149 assert( bufmem != NULL );
152 assert( 0 < i && i <= n );
153 assert( eigenvalue != NULL );
154 assert( ( ! geteigenvectors) || eigenvector != NULL );
157 JOBZ = geteigenvectors ?
'V' :
'N';
172 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
176 &ABSTOL, &M, NULL, NULL,
178 &LWORK, &WISIZE, &LIWORK,
183 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
191 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &WORK, (
int) LWORK) );
192 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &IWORK, (
int) LIWORK) );
193 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &WTMP, (
int) N) );
194 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &ISUPPZ, 2) );
200 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
204 &ABSTOL, &M, WTMP, eigenvector,
206 &LWORK, IWORK, &LIWORK,
211 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
216 *eigenvalue = WTMP[0];
219 BMSfreeBufferMemoryArray(bufmem, &ISUPPZ);
220 BMSfreeBufferMemoryArray(bufmem, &WTMP);
221 BMSfreeBufferMemoryArray(bufmem, &IWORK);
222 BMSfreeBufferMemoryArray(bufmem, &WORK);
233 SCIP_Real* eigenvalues,
234 SCIP_Real* eigenvectors
238 #if ( SDPA_VERSION == 740 )
260 assert( bufmem != NULL );
263 assert( eigenvalues != NULL );
264 assert( eigenvectors != NULL );
280 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
284 &ABSTOL, &M, NULL, NULL,
286 &LWORK, &WISIZE, &LIWORK,
291 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
299 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &WORK, (
int) LWORK) );
300 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &IWORK, (
int) LIWORK) );
301 BMS_CALL( BMSallocBufferMemoryArray(bufmem, &ISUPPZ, (
int) 2 * N) );
307 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
311 &ABSTOL, &M, eigenvalues, eigenvectors,
313 &LWORK, IWORK, &LIWORK,
318 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
323 BMSfreeBufferMemoryArray(bufmem, &ISUPPZ);
324 BMSfreeBufferMemoryArray(bufmem, &IWORK);
325 BMSfreeBufferMemoryArray(bufmem, &WORK);
364 F77_FUNC(dgemv, DGEMV)(&TRANS, &M, &N, &ALPHA, A, &LDA, X, &INCX, &BETA, Y, &INCY);
374 SCIP_Bool transposeA,
378 SCIP_Bool transposeB,
393 assert( (transposeA && transposeB && (nrowsA == ncolsB)) || (transposeA && !transposeB && (nrowsA == nrowsB))
394 || (!transposeA && transposeB && (ncolsA == ncolsB)) || (!transposeA && !transposeB && (ncolsA == nrowsB)) );
396 TRANSA = transposeA ?
'T' :
'N';
397 TRANSB = transposeB ?
'T' :
'N';
398 M = transposeA ? ncolsA : nrowsA;
399 N = transposeB ? nrowsB : ncolsB;
400 K = transposeA ? nrowsA : ncolsA;
402 LDA = transposeA ? K : M;
403 LDB = transposeB ? N : K;
407 F77_FUNC(dgemm, DGEMM)(&TRANSA, &TRANSB, &M, &N, &K, &ALPHA, matrixA, &LDA, matrixB, &LDB, &BETA, result, &LDC);
long long int LAPACKINTTYPE
SCIP_RETCODE SCIPlapackComputeEigenvectorDecomposition(BMS_BUFMEM *bufmem, int n, SCIP_Real *A, SCIP_Real *eigenvalues, SCIP_Real *eigenvectors)
EXTERN SCIP_RETCODE SCIPlapackComputeIthEigenvalue(BMS_BUFMEM *bufmem, SCIP_Bool geteigenvectors, int n, SCIP_Real *A, int i, SCIP_Real *eigenvalue, SCIP_Real *eigenvector)
Macros for calling Fortran-functions.
SCIP_RETCODE SCIPlapackMatrixVectorMult(int nrows, int ncols, SCIP_Real *matrix, SCIP_Real *vector, SCIP_Real *result)
SCIP_RETCODE SCIPlapackMatrixMatrixMult(int nrowsA, int ncolsA, SCIP_Real *matrixA, SCIP_Bool transposeA, int nrowsB, int ncolsB, SCIP_Real *matrixB, SCIP_Bool transposeB, SCIP_Real *result)
#define SCIP_RealTOINT(x)
#define F77_FUNC(name, NAME)
interface methods for eigenvector computation and matrix multiplication using different versions of L...