30 #if defined (INTEL_CXML)
31 #define CHAR_MACRO(char_var) &char_var, 1
33 #define CHAR_MACRO(char_var) &char_var
36 #include "FortranRoutines.h"
40 void FortranRoutines::SCAL_INCX(
int N,
double ALPHA,
double *X,
int incX)
const {
41 DSCAL_F77(&N, &ALPHA, X, &incX);
45 void FortranRoutines::SWAP(
int N,
double *X,
int incx,
double *Y,
int incy)
const {
46 F77_FUNC(dswap,DSWAP)(&N, X, &incx, Y, &incy);
54 void FortranRoutines::GEQRF(
int M,
int N,
double *A,
int lda,
double *tau,
double *work,
55 int lwork,
int *info)
const {
56 F77_FUNC(dgeqrf,DGEQRF)(&M, &N, A, &lda, tau, work, &lwork, info);
60 void FortranRoutines::ORMQR(
char SIDE,
char TRANS,
int M,
int N,
int K,
double *A,
int lda,
61 double *tau,
double *C,
int ldc,
double *work,
int lwork,
63 F77_FUNC(dormqr,DORMQR)(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &M, &N, &K, A, &lda, tau,
64 C, &ldc, work, &lwork, info);
68 void FortranRoutines::SPEV(
char JOBZ,
char UPLO,
int N,
double *A,
double *W,
double *Z,
69 int ldz,
double *work,
int *info)
const {
70 F77_FUNC(dspev,DSPEV)(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &N, A, W, Z, &ldz, work, info);
74 void FortranRoutines::STEQR(
char COMPZ,
int N,
double *D,
double *E,
double *Z,
int ldz,
75 double *work,
int *info)
const {
76 F77_FUNC(dsteqr,DSTEQR)(CHAR_MACRO(COMPZ), &N, D, E, Z, &ldz, work, info);
80 void FortranRoutines::SYEV(
char JOBZ,
char UPLO,
int N,
double *A,
int lda,
double *W,
81 double *work,
int lwork,
int *info)
const {
82 F77_FUNC(dsyev,DSYEV)(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &N, A, &lda, W, work, &lwork, info);
86 void FortranRoutines::SYGV(
int itype,
char JOBZ,
char UPLO,
int N,
double *A,
int lda,
87 double *B,
int ldb,
double *W,
double *work,
int lwork,
89 F77_FUNC(dsygv,DSYGV)(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &N, A, &lda, B, &ldb,
90 W, work, &lwork, info);
94 int FortranRoutines::LAENV(
int ispec,
char *NAME,
char *OPTS,
int N1,
int N2,
int N3,
95 int N4,
int len_name,
int len_opts)
const {
96 #if defined (INTEL_CXML)
97 return F77_FUNC(ilaenv,ILAENV)(&ispec, NAME, len_name, OPTS, len_opts, &N1, &N2, &N3, &N4);
99 return F77_FUNC(ilaenv,ILAENV)(&ispec, NAME, OPTS, &N1, &N2, &N3, &N4, len_name, len_opts);
107 void FortranRoutines::SAUPD(
int *ido,
char BMAT,
int N,
char *which,
int nev,
double tol,
108 double *resid,
int ncv,
double *V,
int ldv,
int *iparam,
109 int *ipntr,
double *workd,
double *workl,
int lworkl,
int *info,
111 #if defined (INTEL_CXML)
112 F77_FUNC(mydsaupd,MYDSAUPD)(ido, &BMAT, 1, &N, which, 2, &nev, &tol, resid, &ncv, V, &ldv,
113 iparam, ipntr, workd, workl, &lworkl, info, &verbose);
115 F77_FUNC(mydsaupd,MYDSAUPD)(ido, &BMAT, &N, which, &nev, &tol, resid, &ncv, V, &ldv,
116 iparam, ipntr, workd, workl, &lworkl, info, &verbose, 1, 2);
121 void FortranRoutines::SEUPD(LOGICAL rvec,
char HOWMNY, LOGICAL *select,
double *D,
122 double *Z,
int ldz,
double sigma,
char BMAT,
int N,
123 char *which,
int nev,
double tol,
double *resid,
int ncv,
double *V,
124 int ldv,
int *iparam,
int *ipntr,
double *workd,
double *workl,
125 int lworkl,
int *info)
const {
126 #if defined (INTEL_CXML)
127 F77_FUNC(dseupd,DSEUPD)(&rvec, &HOWMNY, 1, select, D, Z, &ldz, &sigma, &BMAT, 1, &N,
128 which, 2, &nev, &tol, resid, &ncv, V, &ldv, iparam, ipntr, workd, workl, &lworkl,
131 F77_FUNC(dseupd,DSEUPD)(&rvec, &HOWMNY, select, D, Z, &ldz, &sigma, &BMAT, &N, which,
132 &nev, &tol, resid, &ncv, V, &ldv, iparam, ipntr, workd, workl, &lworkl, info,
142 void FortranRoutines::PSAUPD(MPI_Comm MyComm,
int *ido,
char BMAT,
int N,
char *which,
int nev,
143 double tol,
double *resid,
int ncv,
double *V,
int ldv,
int *iparam,
144 int *ipntr,
double *workd,
double *workl,
int lworkl,
int *info,
146 #if defined (INTEL_CXML)
147 F77_FUNC(mypdsaupd,MYPDSAUPD)(&MyComm, ido, &BMAT, 1, &N, which, 2, &nev, &tol, resid, &ncv,
148 V, &ldv, iparam, ipntr, workd, workl, &lworkl, info, &verbose);
150 F77_FUNC(mypdsaupd,MYPDSAUPD)(&MyComm, ido, &BMAT, &N, which, &nev, &tol, resid, &ncv, V, &ldv,
151 iparam, ipntr, workd, workl, &lworkl, info, &verbose, 1, 2);
156 void FortranRoutines::PSEUPD(MPI_Comm MyComm, LOGICAL rvec,
char HOWMNY, LOGICAL *select,
157 double *D,
double *Z,
int ldz,
double sigma,
char BMAT,
int N,
158 char *which,
int nev,
double tol,
double *resid,
int ncv,
double *V,
159 int ldv,
int *iparam,
int *ipntr,
double *workd,
double *workl,
160 int lworkl,
int *info)
const {
161 #if defined (INTEL_CXML)
162 F77_FUNC(pdseupd,PDSEUPD)(&MyComm, &rvec, &HOWMNY, 1, select, D, Z, &ldz, &sigma, &BMAT, 1, &N,
163 which, 2, &nev, &tol, resid, &ncv, V, &ldv, iparam, ipntr, workd, workl, &lworkl, info);
165 F77_FUNC(pdseupd,PDSEUPD)(&MyComm, &rvec, &HOWMNY, select, D, Z, &ldz, &sigma, &BMAT, &N,
166 which, &nev, &tol, resid, &ncv, V, &ldv, iparam, ipntr, workd, workl, &lworkl, info,