Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
FortranRoutines.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Anasazi: Block Eigensolvers Package
4 //
5 // Copyright 2004 NTESS and the Anasazi contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 // This software is a result of the research described in the report
11 //
12 // "A comparison of algorithms for modal analysis in the absence
13 // of a sparse direct method", P. Arbenz, R. Lehoucq, and U. Hetmaniuk,
14 // Sandia National Laboratories, Technical report SAND2003-1028J.
15 //
16 // It is based on the Epetra, AztecOO, and ML packages defined in the Trilinos
17 // framework ( http://trilinos.org/ ).
18 
19 /* for INTEL_CXML, the second arg may need to be changed to 'one'. If so
20 the appropriate declaration of one will need to be added back into
21 functions that include the macro:
22 #if defined (INTEL_CXML)
23  unsigned int one=1;
24 #endif
25 */
26 
27 #ifdef CHAR_MACRO
28 #undef CHAR_MACRO
29 #endif
30 #if defined (INTEL_CXML)
31 #define CHAR_MACRO(char_var) &char_var, 1
32 #else
33 #define CHAR_MACRO(char_var) &char_var
34 #endif
35 
36 #include "FortranRoutines.h"
37 
38 // Double precision BLAS 1 //
39 
40 void FortranRoutines::SCAL_INCX(int N, double ALPHA, double *X, int incX) const {
41  DSCAL_F77(&N, &ALPHA, X, &incX);
42  return;
43 }
44 
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);
47  return;
48 }
49 
50 
51 // Double precision LAPACK //
52 
53 
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);
57  return;
58 }
59 
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,
62  int *info) const {
63  F77_FUNC(dormqr,DORMQR)(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &M, &N, &K, A, &lda, tau,
64  C, &ldc, work, &lwork, info);
65  return;
66 }
67 
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);
71  return;
72 }
73 
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);
77  return;
78 }
79 
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);
83  return;
84 }
85 
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,
88  int *info) const {
89  F77_FUNC(dsygv,DSYGV)(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &N, A, &lda, B, &ldb,
90  W, work, &lwork, info);
91  return;
92 }
93 
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);
98 #else
99  return F77_FUNC(ilaenv,ILAENV)(&ispec, NAME, OPTS, &N1, &N2, &N3, &N4, len_name, len_opts);
100 #endif
101 }
102 
103 
104 // Double precision ARPACK routines
105 
106 
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,
110  int verbose) const {
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);
114 #else
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);
117 #endif
118  return;
119 }
120 
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,
129  info);
130 #else
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,
133  1, 1, 2);
134 #endif
135  return;
136 }
137 
138 #ifdef EPETRA_MPI
139 
140 // Double precision PARPACK routines
141 
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,
145  int verbose) const {
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);
149 #else
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);
152 #endif
153  return;
154 }
155 
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);
164 #else
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,
167  1, 1, 2);
168 #endif
169  return;
170 }
171 
172 #endif
173