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 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: #pragma once

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

 30: /* SUBMANSEC = Sys */

 32: /*MC
 33:    PetscCallBLAS - Calls a BLAS or LAPACK routine so that the stack trace returned from any signal received inside the function call
 34:    includes the name of the BLAS/LAPACK routine

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

 40:    Not Collective

 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 Notes:
 49:    This does not check error codes returned from the BLAS/LAPACK routine or ever return from the current subroutine. It merely pushes onto the PETSc
 50:    stack the name of the BLAS/LAPACK routine before calling the routine and removes it after a successful call.

 52:    LAPACK routines that return info codes should be followed by
 53: .vb
 54:    PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, ...)
 55: .ve

 57:    This macro exists so that when a BLAS/LAPACK routine results in a crash or corrupts memory, they get blamed instead of PETSc.

 59: .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscStackCallExternalVoid()`
 60: M*/
 61: #define PetscCallBLAS(name, routine) \
 62:   do { \
 63:     PetscStackPushExternal(name); \
 64:     routine; \
 65:     PetscStackPop; \
 66:   } while (0)

 68: static inline void PetscMissingLapack(const char *fname, ...)
 69: {
 70:   SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_SUP, "%s - Lapack routine is unavailable.", fname);
 71: }

 73: #include <petscblaslapack_mangle.h>

 75: BLAS_EXTERN void LAPACKgetrf_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
 76: BLAS_EXTERN void LAPACKREALgetrf_(const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *, PetscBLASInt *);
 77: BLAS_EXTERN void LAPACKgetri_(const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 78: BLAS_EXTERN void LAPACKREALgetri_(const PetscBLASInt *, PetscReal *, const PetscBLASInt *, const PetscBLASInt *, PetscReal *, const PetscBLASInt *, PetscBLASInt *);
 79: #if !defined(PETSC_MISSING_LAPACK_ORGQR)
 80: BLAS_EXTERN void LAPACKorgqr_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 81: #else
 82:   #define LAPACKorgqr_(a, b, c, d, e, f, g, h, i) PetscMissingLapack("ORGQR", a, b, c, d, e, f, g, h, i)
 83: #endif
 84: BLAS_EXTERN void LAPACKgeqrf_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
 85: #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SNRM2_RETURNS_DOUBLE)
 86: BLAS_EXTERN double BLASnrm2_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 87: #else
 88: BLAS_EXTERN PetscReal BLASnrm2_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 89: #endif
 90: BLAS_EXTERN void BLASscal_(const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
 91: BLAS_EXTERN void BLAScopy_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 92: BLAS_EXTERN void BLASswap_(const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 93: BLAS_EXTERN void BLASaxpy_(const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
 94: #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SNRM2_RETURNS_DOUBLE)
 95: BLAS_EXTERN double BLASasum_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 96: #else
 97: BLAS_EXTERN PetscReal BLASasum_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
 98: #endif
 99: BLAS_EXTERN void LAPACKpttrf_(const PetscBLASInt *, PetscReal *, PetscScalar *, PetscBLASInt *);
100: #if !defined(PETSC_MISSING_LAPACK_STEIN)
101: 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 *);
102: #else
103:   #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)
104: #endif
105: BLAS_EXTERN void LAPACKgesv_(const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);

107: BLAS_EXTERN void LAPACKtrtri_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
108: BLAS_EXTERN void LAPACKpotrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
109: BLAS_EXTERN void LAPACKpotri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
110: BLAS_EXTERN void LAPACKpotrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
111: BLAS_EXTERN void LAPACKsytrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
112: BLAS_EXTERN void LAPACKsytrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
113: #if !defined(PETSC_MISSING_LAPACK_SYTRI)
114: BLAS_EXTERN void LAPACKsytri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, PetscBLASInt *);
115: #else
116:   #define LAPACKsytri_(a, b, c, d, e, f, g) PetscMissingLapack("SYTRI", a, b, c, d, e, f, g)
117: #endif
118: BLAS_EXTERN void BLASsyrk_(const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, PetscScalar *, const PetscBLASInt *);
119: 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 *);
120: 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 *);
121: BLAS_EXTERN void LAPACKgetrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
122: BLAS_EXTERN void BLAStrmv_(const char *, const char *, const char *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
123: BLAS_EXTERN void BLAStrsv_(const char *, const char *, const char *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *);
124: 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 *);
125: 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 *);
126: 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 *);
127: 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 *);
128: #if !defined(PETSC_MISSING_LAPACK_ORMQR)
129: BLAS_EXTERN void LAPACKormqr_(const char *, const char *, PetscBLASInt *, PetscBLASInt *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscScalar *, PetscBLASInt *, PetscScalar *, PetscBLASInt *, PetscBLASInt *);
130: #else
131:   #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)
132: #endif
133: #if !defined(PETSC_MISSING_LAPACK_STEGR)
134: 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 *);
135: #else
136:   #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)
137: #endif
138: #if !defined(PETSC_MISSING_LAPACK_STEQR)
139: BLAS_EXTERN void LAPACKsteqr_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
140: BLAS_EXTERN void LAPACKREALsteqr_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
141: #else
142:   #define LAPACKsteqr_(a, b, c, d, e, f, g, h)     PetscMissingLapack("STEQR", a, b, c, d, e, f, g, h)
143:   #define LAPACKREALsteqr_(a, b, c, d, e, f, g, h) PetscMissingLapack("STEQR", a, b, c, d, e, f, g, h)
144: #endif
145: #if !defined(PETSC_MISSING_LAPACK_STEV)
146: BLAS_EXTERN void LAPACKstev_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
147: BLAS_EXTERN void LAPACKREALstev_(const char *, const PetscBLASInt *, PetscReal *, PetscReal *, PetscReal *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
148: #else
149:   #define LAPACKstev_(a, b, c, d, e, f, g, h)     PetscMissingLapack("STEV", a, b, c, d, e, f, g, h)
150:   #define LAPACKREALstev_(a, b, c, d, e, f, g, h) PetscMissingLapack("STEV", a, b, c, d, e, f, g, h)
151: #endif
152: #if !defined(PETSC_MISSING_LAPACK_HGEQZ)
153: 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 *);
154: #else
155:   #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)
156: #endif
157: #if !defined(PETSC_MISSING_LAPACK_TRTRS)
158: BLAS_EXTERN void LAPACKtrtrs_(const char *, const char *, const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
159: #else
160:   #define LAPACKtrtrs_(a, b, c, d, e, f, g, h, i, j) PetscMissingLapack("TRTRS", a, b, c, d, e, f, g, h, i, j)
161: #endif
162: BLAS_EXTERN void LAPACKgels_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);

164: /* handle complex dot() with special code */
165: #if defined(PETSC_USE_COMPLEX)
166: static inline PetscScalar BLASdot_(const PetscBLASInt *n, const PetscScalar *x, const PetscBLASInt *sx, const PetscScalar *y, const PetscBLASInt *sy)
167: {
168:   PetscScalar sum = 0.0;
169:   PetscInt    i, j, k;
170:   if (*sx == 1 && *sy == 1) {
171:     for (i = 0; i < *n; i++) sum += PetscConj(x[i]) * y[i];
172:   } else {
173:     for (i = 0, j = 0, k = 0; i < *n; i++, j += *sx, k += *sy) sum += PetscConj(x[j]) * y[k];
174:   }
175:   return sum;
176: }
177: static inline PetscScalar BLASdotu_(const PetscBLASInt *n, const PetscScalar *x, const PetscBLASInt *sx, const PetscScalar *y, const PetscBLASInt *sy)
178: {
179:   PetscScalar sum = 0.0;
180:   PetscInt    i, j, k;
181:   if (*sx == 1 && *sy == 1) {
182:     for (i = 0; i < *n; i++) sum += x[i] * y[i];
183:   } else {
184:     for (i = 0, j = 0, k = 0; i < *n; i++, j += *sx, k += *sy) sum += x[j] * y[k];
185:   }
186:   return sum;
187: }
188: #else
189:   #if defined(PETSC_USE_REAL_SINGLE) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
190: BLAS_EXTERN double BLASdot_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
191: BLAS_EXTERN double BLASdotu_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
192:   #else
193: BLAS_EXTERN PetscScalar BLASdot_(const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *);
194:   #endif
195: #endif

197: /* Some functions prototypes do not exist for reals */
198: #if defined(PETSC_USE_COMPLEX)
199: BLAS_EXTERN void LAPACKhetrf_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
200: BLAS_EXTERN void LAPACKhetrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
201: BLAS_EXTERN void LAPACKhetri_(const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, PetscBLASInt *);
202: BLAS_EXTERN void LAPACKheev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
203: #endif
204: /* Some functions prototypes differ between real and complex */
205: #if defined(PETSC_USE_COMPLEX)
206:   #if !defined(PETSC_MISSING_LAPACK_GELSS)
207: 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 *);
208:   #else
209:     #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)
210:   #endif
211: BLAS_EXTERN void LAPACKsyev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
212: 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 *);
213: BLAS_EXTERN void LAPACKsygv_(const PetscBLASInt *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscBLASInt *);
214: 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 *);
215: BLAS_EXTERN void LAPACKpttrs_(const char *, const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
216:   #if !defined(PETSC_MISSING_LAPACK_GERFS)
217: 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 *);
218:   #else
219:     #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)
220:   #endif
221:   #if !defined(PETSC_MISSING_LAPACK_TRSEN)
222: 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 *);
223:   #else
224:     #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)
225:   #endif
226:   #if !defined(PETSC_MISSING_LAPACK_TGSEN)
227: 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 *);
228:   #else
229:     #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)
230:   #endif
231:   #if !defined(PETSC_MISSING_LAPACK_GGES)
232: 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 *);
233:   #else
234:     #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)
235:   #endif
236:   #if !defined(PETSC_MISSING_LAPACK_HSEQR)
237: 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 *);
238:   #else
239:     #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)
240:   #endif
241: #else /* !defined(PETSC_USE_COMPLEX) */
242:   #if !defined(PETSC_MISSING_LAPACK_GELSS)
243: BLAS_EXTERN void LAPACKgelss_(const PetscBLASInt *, const PetscBLASInt *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, const PetscReal *, PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
244:   #else
245:     #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)
246:   #endif
247: BLAS_EXTERN void LAPACKsyev_(const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscReal *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
248: 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 *);
249: BLAS_EXTERN void LAPACKsygv_(const PetscBLASInt *, const char *, const char *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, const PetscBLASInt *, PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
250: 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 *);
251: BLAS_EXTERN void LAPACKpttrs_(const PetscBLASInt *, const PetscBLASInt *, const PetscReal *, const PetscScalar *, PetscScalar *, const PetscBLASInt *, PetscBLASInt *);
252:   #if !defined(PETSC_MISSING_LAPACK_STEBZ)
253: 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 *);
254:   #else
255:     #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)
256:   #endif
257:   #if !defined(PETSC_MISSING_LAPACK_GERFS)
258: 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 *);
259:   #else
260:     #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)
261:   #endif
262:   #if !defined(PETSC_MISSING_LAPACK_TRSEN)
263: 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 *);
264:   #else
265:     #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)
266:   #endif
267:   #if !defined(PETSC_MISSING_LAPACK_TGSEN)
268: 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 *);
269:   #else
270:     #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)
271:   #endif
272:   #if !defined(PETSC_MISSING_LAPACK_GGES)
273: 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 *);
274:   #else
275:     #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)
276:   #endif
277:   #if !defined(PETSC_MISSING_LAPACK_HSEQR)
278: 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 *);
279:   #else
280:     #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)
281:   #endif
282: #endif /* defined(PETSC_USE_COMPLEX) */

284: #if defined(PETSC_USE_COMPLEX)
285: 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 *);
286: 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 *);
287: #else
288: 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 *);
289: 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 *);
290: #endif