Actual source code: petscblaslapack.h

  1: /*
  2:   This file dispatches between various header files for blas/lapack distributions to handle the name mangling.
  3:   It also provides C prototypes for all the BLAS/LAPACK functions that PETSc uses

  5:   This is not included automatically by petscsys.h because some external packages include their own prototypes for
  6:   certain BLAS/LAPACK functions that conflict with the ones given here. Hence this should only be included when needed.

  8:   The BLAS/LAPACK name mangling is almost (but not always) the same as the Fortran mangling; and exists even if there is
  9:   not Fortran compiler.

 11:   PETSC_BLASLAPACK_UNDERSCORE BLAS/LAPACK function have an underscore at the end of each function name
 12:   PETSC_BLASLAPACK_CAPS BLAS/LAPACK function names are all in capital letters
 13:   PETSC_BLASLAPACK_C BLAS/LAPACK function names have no mangling

 15:   PETSC_BLASLAPACK_SINGLEISDOUBLE - for Cray systems where the BLAS/LAPACK single precision (i.e. Fortran single precision is actually 64 bits)
 16:                                     old Cray vector machines used to be this way, it is is not clear if any exist now.

 18:   PetscBLASInt is almost always 32 bit integers but can be 64 bit integers for certain usages of MKL and OpenBLAS BLAS/LAPACK libraries

 20: */
 21: #ifndef _BLASLAPACK_H
 22: #define _BLASLAPACK_H

 24: #include <petscconf.h>
 25: #if defined(__cplusplus)
 26:   #define BLAS_EXTERN extern "C"
 27: #else
 28:   #define BLAS_EXTERN extern
 29: #endif

 31: /* SUBMANSEC = Sys */

 33: /*MC
 34:     PetscCallBLAS - Calls a BLAS or LAPACK routine with error check handling

 36:     Not collective

 38:     Synopsis:
 39: #include <petscsys.h>
 40:    void PetscCallBLAS(char *name,routine)

 42:   Input Parameters:
 43: +   name - string that gives the name of the function being called
 44: -   routine - actual call to the routine including its arguments

 46:    Level: developer

 48:    Developer Note:
 49:    This is so that when a user or external library routine results in a crash or corrupts memory, they get blamed instead of PETSc.

 51: .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscStackCallExternalVoid()`
 52: M*/
 53: #define PetscCallBLAS(name, routine) \
 54:   do { \
 55:     PetscStackPushExternal(name); \
 56:     routine; \
 57:     PetscStackPop; \
 58:   } while (0)

 60: static inline void PetscMissingLapack(const char *fname, ...)
 61: {
 62:   PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL, "%s - Lapack routine is unavailable.", fname);
 63:   MPI_Abort(PETSC_COMM_SELF, PETSC_ERR_SUP);
 64: }

 66: #include <petscblaslapack_mangle.h>

 68: BLAS_EXTERN void LAPACKgetrf_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
 69: BLAS_EXTERN void LAPACKREALgetrf_(const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
 70: BLAS_EXTERN void LAPACKgetri_(const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 71: BLAS_EXTERN void LAPACKREALgetri_(const PetscBLASInt *, PetscReal *, const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *);
 72: #if !defined(PETSC_MISSING_LAPACK_ORGQR)
 73: BLAS_EXTERN void LAPACKorgqr_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 74: #else
 75:   #define LAPACKorgqr_(a, b, c, d, e, f, g, h, i) PetscMissingLapack("ORGQR", a, b, c, d, e, f, g, h, i)
 76: #endif
 77: BLAS_EXTERN void LAPACKgeqrf_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 78: #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SNRM2_RETURNS_DOUBLE)
 79: BLAS_EXTERN double BLASnrm2_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 80: #else
 81: BLAS_EXTERN PetscReal BLASnrm2_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 82: #endif
 83: BLAS_EXTERN void BLASscal_(const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
 84: BLAS_EXTERN void BLAScopy_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 85: BLAS_EXTERN void BLASswap_(const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 86: BLAS_EXTERN void BLASaxpy_(const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 87: #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SNRM2_RETURNS_DOUBLE)
 88: BLAS_EXTERN double BLASasum_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 89: #else
 90: BLAS_EXTERN PetscReal BLASasum_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 91: #endif
 92: BLAS_EXTERN void LAPACKpttrf_(const PetscBLASInt *, PetscReal *, PetscScalar *, PetscBLASInt *);
 93: #if !defined(PETSC_MISSING_LAPACK_STEIN)
 94: BLAS_EXTERN void LAPACKstein_(const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
 95: #else
 96:   #define LAPACKstein_(a, b, c, d, e, f, g, h, i, j, k, l, m) PetscMissingLapack("STEIN", a, b, c, d, e, f, g, h, i, j, k, l)
 97: #endif
 98: BLAS_EXTERN void LAPACKgesv_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);

100: BLAS_EXTERN void LAPACKpotrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
101: BLAS_EXTERN void LAPACKpotri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
102: BLAS_EXTERN void LAPACKpotrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
103: BLAS_EXTERN void LAPACKsytrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
104: BLAS_EXTERN void LAPACKsytrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
105: #if !defined(PETSC_MISSING_LAPACK_SYTRI)
106: BLAS_EXTERN void LAPACKsytri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, PetscBLASInt *);
107: #else
108:   #define LAPACKsytri_(a, b, c, d, e, f, g) PetscMissingLapack("SYTRI", a, b, c, d, e, f, g)
109: #endif
110: BLAS_EXTERN void BLASsyrk_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
111: BLAS_EXTERN void BLASsyr2k_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
112: BLAS_EXTERN void BLASgemv_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
113: BLAS_EXTERN void LAPACKgetrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
114: BLAS_EXTERN void BLAStrmv_(const char *, const char *, const char *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
115: BLAS_EXTERN void BLASgemm_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
116: BLAS_EXTERN void BLASREALgemm_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscReal *, const PetscBLASInt *, const PetscReal *, PetscReal *, const PetscBLASInt *);
117: BLAS_EXTERN void BLASsymm_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
118: BLAS_EXTERN void BLAStrsm_(const char *, const char *, const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
119: #if !defined(PETSC_MISSING_LAPACK_ORMQR)
120: BLAS_EXTERN void LAPACKormqr_(const char *, const char *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscBLASInt *);
121: #else
122:   #define LAPACKormqr_(a, b, c, d, e, f, g, h, i, j, k, l, m) PetscMissingLapack("ORMQR", a, b, c, d, e, f, g, h, i, j, k, l, m)
123: #endif
124: #if !defined(PETSC_MISSING_LAPACK_STEGR)
125: BLAS_EXTERN void LAPACKstegr_(const char *, const char *, const PetscBLASInt *, PetscReal *, PetscReal *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, const PetscBLASInt *, PetscBLASInt *);
126: #else
127:   #define LAPACKstegr_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) PetscMissingLapack("STEGR", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t)
128: #endif
129: #if !defined(PETSC_MISSING_LAPACK_STEQR)
130: BLAS_EXTERN void LAPACKsteqr_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
131: BLAS_EXTERN void LAPACKREALsteqr_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
132: #else
133:   #define LAPACKsteqr_(a, b, c, d, e, f, g, h)     PetscMissingLapack("STEQR", a, b, c, d, e, f, g, h)
134:   #define LAPACKREALsteqr_(a, b, c, d, e, f, g, h) PetscMissingLapack("STEQR", a, b, c, d, e, f, g, h)
135: #endif
136: #if !defined(PETSC_MISSING_LAPACK_HGEQZ)
137: BLAS_EXTERN void LAPACKhgeqz_(const char *, const char *, const char *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscBLASInt *);
138: #else
139:   #define LAPACKhgeqz_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) PetscMissingLapack("HGEQZ", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t)
140: #endif
141: #if !defined(PETSC_MISSING_LAPACK_TRTRS)
142: BLAS_EXTERN void LAPACKtrtrs_(const char *, const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
143: #else
144:   #define LAPACKtrtrs_(a, b, c, d, e, f, g, h, i, j) PetscMissingLapack("TRTRS", a, b, c, d, e, f, g, h, i, j)
145: #endif
146: BLAS_EXTERN void LAPACKgels_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);

148: /* handle complex dot() with special code */
149: #if defined(PETSC_USE_COMPLEX)
150: static inline PetscScalar BLASdot_(const PetscBLASInt *n, const PetscScalar *x, const PetscBLASInt *sx, const PetscScalar *y, const PetscBLASInt *sy)
151: {
152:   PetscScalar sum = 0.0;
153:   PetscInt    i, j, k;
154:   if (*sx == 1 && *sy == 1) {
155:     for (i = 0; i < *n; i++) sum += PetscConj(x[i]) * y[i];
156:   } else {
157:     for (i = 0, j = 0, k = 0; i < *n; i++, j += *sx, k += *sy) sum += PetscConj(x[j]) * y[k];
158:   }
159:   return sum;
160: }
161: static inline PetscScalar BLASdotu_(const PetscBLASInt *n, const PetscScalar *x, const PetscBLASInt *sx, const PetscScalar *y, const PetscBLASInt *sy)
162: {
163:   PetscScalar sum = 0.0;
164:   PetscInt    i, j, k;
165:   if (*sx == 1 && *sy == 1) {
166:     for (i = 0; i < *n; i++) sum += x[i] * y[i];
167:   } else {
168:     for (i = 0, j = 0, k = 0; i < *n; i++, j += *sx, k += *sy) sum += x[j] * y[k];
169:   }
170:   return sum;
171: }
172: #else
173:   #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
174: BLAS_EXTERN double    BLASdot_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
175: BLAS_EXTERN double    BLASdotu_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
176:   #else
177: BLAS_EXTERN PetscScalar BLASdot_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
178:   #endif
179: #endif

181: /* Some functions prototypes do not exist for reals */
182: #if defined(PETSC_USE_COMPLEX)
183: BLAS_EXTERN void LAPACKhetrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
184: BLAS_EXTERN void LAPACKhetrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
185: BLAS_EXTERN void LAPACKhetri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, PetscBLASInt *);
186: BLAS_EXTERN void LAPACKheev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
187: #endif
188: /* Some functions prototypes differ between real and complex */
189: #if defined(PETSC_USE_COMPLEX)
190:   #if !defined(PETSC_MISSING_LAPACK_GELSS)
191: BLAS_EXTERN void LAPACKgelss_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, const PetscReal *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
192:   #else
193:     #define LAPACKgelss_(a, b, c, d, e, f, g, h, i, j, k, l, m, n) PetscMissingLapack("GELSS", a, b, c, d, e, f, g, h, i, j, k, l, m, n)
194:   #endif
195: BLAS_EXTERN void LAPACKsyev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
196: BLAS_EXTERN void LAPACKsyevx_(const char *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
197: BLAS_EXTERN void LAPACKsygv_(const PetscBLASInt *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
198: BLAS_EXTERN void LAPACKsygvx_(PetscBLASInt *, const char *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
199: BLAS_EXTERN void LAPACKpttrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
200:   #if !defined(PETSC_MISSING_LAPACK_GERFS)
201: BLAS_EXTERN void LAPACKgerfs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, PetscReal *, PetscBLASInt *);
202:   #else
203:     #define LAPACKgerfs_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q) PetscMissingLapack("GERFS", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q)
204:   #endif
205:   #if !defined(PETSC_MISSING_LAPACK_TRSEN)
206: BLAS_EXTERN void LAPACKtrsen_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
207:   #else
208:     #define LAPACKtrsen_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o) PetscMissingLapack("TRSEN", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o)
209:   #endif
210:   #if !defined(PETSC_MISSING_LAPACK_TGSEN)
211: BLAS_EXTERN void LAPACKtgsen_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, const PetscBLASInt *, PetscBLASInt *);
212:   #else
213:     #define LAPACKtgsen_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u, v, w, x) PetscMissingLapack("TGSEN", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u, v, w, x)
214:   #endif
215:   #if !defined(PETSC_MISSING_LAPACK_GGES)
216: BLAS_EXTERN void LAPACKgges_(const char *, const char *, const char *, PetscBLASInt (*)(void), const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *);
217:   #else
218:     #define LAPACKgges_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u) PetscMissingLapack("GGES", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u)
219:   #endif
220:   #if !defined(PETSC_MISSING_LAPACK_HSEQR)
221: BLAS_EXTERN void LAPACKhseqr_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
222:   #else
223:     #define LAPACKhseqr_(a, b, c, d, e, f, g, h, i, j, k, l, m) PetscMissingLapack("HSEQR", a, b, c, d, e, f, g, h, i, j, k, l, m)
224:   #endif
225: #else /* !defined(PETSC_USE_COMPLEX) */
226:   #if !defined(PETSC_MISSING_LAPACK_GELSS)
227: BLAS_EXTERN void LAPACKgelss_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, const PetscReal *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
228:   #else
229:     #define LAPACKgelss_(a, b, c, d, e, f, g, h, i, j, k, l, m) PetscMissingLapack("GELSS", a, b, c, d, e, f, g, h, i, j, k, l, m)
230:   #endif
231: BLAS_EXTERN void LAPACKsyev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
232: BLAS_EXTERN void LAPACKsyevx_(const char *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
233: BLAS_EXTERN void LAPACKsygv_(const PetscBLASInt *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
234: BLAS_EXTERN void LAPACKsygvx_(const PetscBLASInt *, const char *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
235: BLAS_EXTERN void LAPACKpttrs_(const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
236:   #if !defined(PETSC_MISSING_LAPACK_STEBZ)
237: BLAS_EXTERN void LAPACKstebz_(const char *, const char *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscReal *, const PetscReal *, PetscBLASInt *, PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *, PetscReal *, PetscBLASInt *, PetscBLASInt *);
238:   #else
239:     #define LAPACKstebz_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r) PetscMissingLapack("STEBZ", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r)
240:   #endif
241:   #if !defined(PETSC_MISSING_LAPACK_GERFS)
242: BLAS_EXTERN void LAPACKgerfs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, PetscBLASInt *, PetscBLASInt *);
243:   #else
244:     #define LAPACKgerfs_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q) PetscMissingLapack("GERFS", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q)
245:   #endif
246:   #if !defined(PETSC_MISSING_LAPACK_TRSEN)
247: BLAS_EXTERN void LAPACKtrsen_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, const PetscBLASInt *, PetscBLASInt *);
248:   #else
249:     #define LAPACKtrsen_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r) PetscMissingLapack("TRSEN", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r)
250:   #endif
251:   #if !defined(PETSC_MISSING_LAPACK_TGSEN)
252: BLAS_EXTERN void LAPACKtgsen_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, const PetscBLASInt *, PetscBLASInt *);
253:   #else
254:     #define LAPACKtgsen_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u, v, w, x, y) PetscMissingLapack("TGSEN", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u, v, w, x, y)
255:   #endif
256:   #if !defined(PETSC_MISSING_LAPACK_GGES)
257: BLAS_EXTERN void LAPACKgges_(const char *, const char *, const char *, PetscBLASInt (*)(void), const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
258:   #else
259:     #define LAPACKgges_(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u) PetscMissingLapack("GGES", a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t, u)
260:   #endif
261:   #if !defined(PETSC_MISSING_LAPACK_HSEQR)
262: BLAS_EXTERN void LAPACKhseqr_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
263:   #else
264:     #define LAPACKhseqr_(a, b, c, d, e, f, g, h, i, j, k, l, m, n) PetscMissingLapack("HSEQR", a, b, c, d, e, f, g, h, i, j, k, l, m, n)
265:   #endif
266: #endif /* defined(PETSC_USE_COMPLEX) */

268: #if defined(PETSC_USE_COMPLEX)
269: BLAS_EXTERN void LAPACKgeev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
270: BLAS_EXTERN void LAPACKgesvd_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
271: #else
272: BLAS_EXTERN void LAPACKgeev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
273: BLAS_EXTERN void LAPACKgesvd_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
274: #endif

276: #endif