49 #include "blockmemshell/memory.h"
50 #include "scip/type_retcode.h"
56 #define BMS_CALL(x) do \
60 SCIPerrorMessage("No memory in function call\n"); \
61 return SCIP_NOMEMORY; \
67 #define SCIP_RealTOINT(x) ((int) (x + 0.5))
78 char* JOBZ,
char* RANGE,
char* UPLO,
79 int* N, SCIP_Real* A,
int* LDA,
80 SCIP_Real* VL, SCIP_Real* VU,
82 SCIP_Real* ABSTOL,
int* M, SCIP_Real* W, SCIP_Real* Z,
83 int* LDZ,
int* ISUPPZ, SCIP_Real* WORK,
84 int* LWORK,
int* IWORK,
int* LIWORK,
89 void F77_FUNC(dgemv, DGEMV)(
char* TRANS,
int* M,
int* N, SCIP_Real* ALPHA, SCIP_Real* A,
int* LDA, SCIP_Real* X,
int* INCX, SCIP_Real* BETA, SCIP_Real* Y,
int* INCY);
105 SCIP_Bool geteigenvectors,
109 SCIP_Real* eigenvalue,
110 SCIP_Real* eigenvector
136 assert( blkmem != NULL );
139 assert( 0 < i && i <= n );
140 assert( eigenvalue != NULL );
141 assert( ( ! geteigenvectors) || eigenvector != NULL );
144 JOBZ = geteigenvectors ?
'V' :
'N';
159 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
163 &ABSTOL, &M, NULL, NULL,
165 &LWORK, &WISIZE, &LIWORK,
170 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
178 BMS_CALL( BMSallocBlockMemoryArray(blkmem, &WORK, LWORK) );
179 BMS_CALL( BMSallocBlockMemoryArray(blkmem, &IWORK, LIWORK) );
180 BMS_CALL( BMSallocBlockMemoryArray(blkmem, &WTMP, N) );
181 BMS_CALL( BMSallocBlockMemoryArray(blkmem, &ISUPPZ, 2) );
187 F77_FUNC(dsyevr, DSYEVR)( &JOBZ, &RANGE, &UPLO,
191 &ABSTOL, &M, WTMP, eigenvector,
193 &LWORK, IWORK, &LIWORK,
198 SCIPerrorMessage(
"There was an error when calling DSYEVR. INFO = %d\n", INFO);
203 *eigenvalue = WTMP[0];
206 BMSfreeBlockMemoryArray(blkmem, &ISUPPZ, 2);
207 BMSfreeBlockMemoryArray(blkmem, &WTMP, N);
208 BMSfreeBlockMemoryArray(blkmem, &IWORK, LIWORK);
209 BMSfreeBlockMemoryArray(blkmem, &WORK, LWORK);
247 F77_FUNC(dgemv, DGEMV)(&TRANS, &M, &N, &ALPHA, A, &LDA, X, &INCX, &BETA, Y, &INCY);
Macros for calling Fortran-functions.
SCIP_RETCODE SCIPlapackMatrixVectorMult(int nrows, int ncols, SCIP_Real *matrix, SCIP_Real *vector, SCIP_Real *result)
#define F77_FUNC(name, NAME)
interface methods for eigenvector computation and matrix multiplication using different versions of L...
#define SCIP_RealTOINT(x)
EXTERN SCIP_RETCODE SCIPlapackComputeIthEigenvalue(BMS_BLKMEM *blkmem, SCIP_Bool geteigenvectors, int n, SCIP_Real *A, int i, SCIP_Real *eigenvalue, SCIP_Real *eigenvector)