42 #ifndef _TEUCHOS_LAPACK_HPP_
43 #define _TEUCHOS_LAPACK_HPP_
89 struct UndefinedLAPACKRoutine
92 static inline T notDefined() {
return T::LAPACK_routine_not_defined_for_this_type(); }
95 template<
typename OrdinalType,
typename ScalarType>
119 void PTTRF(
const OrdinalType& n, ScalarType* d, ScalarType* e, OrdinalType* info)
const;
122 void PTTRS(
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* d,
const ScalarType* e, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
125 void POTRF(
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* info)
const;
128 void POTRS(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
131 void POTRI(
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* info)
const;
135 void POCON(
const char& UPLO,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
138 void POSV(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
141 void POEQU(
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda, MagnitudeType* S, MagnitudeType* scond, MagnitudeType* amax, OrdinalType* info)
const;
144 void PORFS(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
148 TEUCHOS_DEPRECATED
void POSVX(
const char& FACT,
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf,
const char& EQUED, ScalarType* S, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
149 void POSVX(
const char& FACT,
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf,
char* EQUED, ScalarType* S, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
156 void GELS(
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
191 void GELSS(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* S,
const MagnitudeType rcond, OrdinalType* rank, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const;
194 void GELSS(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* S,
const ScalarType& rcond, OrdinalType* rank, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
197 void GGLSE(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& p, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* C, ScalarType* D, ScalarType* X, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
200 void GEQRF (
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
203 void GEQR2 (
const OrdinalType& m,
const OrdinalType& n, ScalarType A[],
const OrdinalType& lda, ScalarType TAU[], ScalarType WORK[], OrdinalType*
const info)
const;
206 void GETRF(
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV, OrdinalType* info)
const;
209 void GETRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
212 void LASCL(
const char& TYPE,
const OrdinalType& kl,
const OrdinalType& ku,
const MagnitudeType cfrom,
const MagnitudeType cto,
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* info)
const;
216 GEQP3(
const OrdinalType& m,
217 const OrdinalType& n, ScalarType* A,
218 const OrdinalType& lda,
222 const OrdinalType& lwork,
223 MagnitudeType* RWORK,
224 OrdinalType* info )
const;
228 LASWP (
const OrdinalType& N,
230 const OrdinalType& LDA,
231 const OrdinalType& K1,
232 const OrdinalType& K2,
233 const OrdinalType IPIV[],
234 const OrdinalType& INCX)
const;
237 void GBTRF(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku, ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV, OrdinalType* info)
const;
240 void GBTRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
243 void GTTRF(
const OrdinalType& n, ScalarType* dl, ScalarType* d, ScalarType* du, ScalarType* du2, OrdinalType* IPIV, OrdinalType* info)
const;
246 void GTTRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* dl,
const ScalarType* d,
const ScalarType* du,
const ScalarType* du2,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
249 void GETRI(
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
256 LATRS (
const char& UPLO,
260 const OrdinalType& N,
262 const OrdinalType& LDA,
264 MagnitudeType* SCALE,
265 MagnitudeType* CNORM,
266 OrdinalType* INFO)
const;
269 void GECON(
const char& NORM,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
272 void GBCON(
const char& NORM,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
278 void GESV(
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
281 void GEEQU(
const OrdinalType& m,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda, ScalarType* R, ScalarType* C, ScalarType* rowcond, ScalarType* colcond, ScalarType* amax, OrdinalType* info)
const;
284 void GERFS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const OrdinalType* IPIV,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
287 void GBEQU(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const ScalarType* A,
const OrdinalType& lda, MagnitudeType* R, MagnitudeType* C, MagnitudeType* rowcond, MagnitudeType* colcond, MagnitudeType* amax, OrdinalType* info)
const;
290 void GBRFS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const OrdinalType* IPIV,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
295 TEUCHOS_DEPRECATED
void GESVX(
const char& FACT,
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf, OrdinalType* IPIV,
const char& EQUED, ScalarType* R, ScalarType* C, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
296 void GESVX(
const char& FACT,
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf, OrdinalType* IPIV,
char* EQUED, ScalarType* R, ScalarType* C, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const;
301 void SYTRD(
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* D, ScalarType* E, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
304 void GEHRD(
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* A,
const OrdinalType& lda, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
307 void TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const;
310 void TRTRI(
const char& UPLO,
const char& DIAG,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda, OrdinalType* info)
const;
318 void SPEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* AP, ScalarType* W, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK, OrdinalType* info)
const;
323 void SYEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* W, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
328 void SYGV(
const OrdinalType& itype,
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* W, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
333 void HEEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* W, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const;
338 void HEGV(
const OrdinalType& itype,
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* W, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType *RWORK, OrdinalType* info)
const;
341 void STEQR(
const char& COMPZ,
const OrdinalType& n, ScalarType* D, ScalarType* E, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK, OrdinalType* info)
const;
346 void HSEQR(
const char& JOB,
const char& COMPZ,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* H,
const OrdinalType& ldh, ScalarType* WR, ScalarType* WI, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
352 void GEES(
const char& JOBVS,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType*, ScalarType*),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, ScalarType* WR, ScalarType* WI, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* BWORK, OrdinalType* info)
const;
357 void GEES(
const char& JOBVS,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType*),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, ScalarType* W, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* BWORK, OrdinalType* info)
const;
362 void GEES(
const char& JOBVS,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, MagnitudeType* WR, MagnitudeType* WI, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* BWORK, OrdinalType* info)
const;
369 void GEEV(
const char& JOBVL,
const char& JOBVR,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* WR, MagnitudeType* WI, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const;
375 void GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* WR, ScalarType* WI, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* SCALE, MagnitudeType* abnrm, MagnitudeType* RCONDE, MagnitudeType* RCONDV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK, OrdinalType* info)
const;
381 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* ALPHAR, MagnitudeType* ALPHAI, ScalarType* BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* lscale, MagnitudeType* rscale, MagnitudeType* abnrm, MagnitudeType* bbnrm, MagnitudeType* RCONDE, MagnitudeType* RCONDV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK, OrdinalType* BWORK, OrdinalType* info)
const;
386 void GGEV(
const char& JOBVL,
const char& JOBVR,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType* BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
392 void TRSEN(
const char& JOB,
const char& COMPQ,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* T,
const OrdinalType& ldt, ScalarType* Q,
const OrdinalType& ldq, MagnitudeType *WR, MagnitudeType *WI, OrdinalType* M, ScalarType* S, MagnitudeType *SEP, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK,
const OrdinalType& liwork, OrdinalType* info )
const;
398 void TGSEN(
const OrdinalType& ijob,
const OrdinalType& wantq,
const OrdinalType& wantz,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType* Q,
const OrdinalType& ldq, ScalarType* Z,
const OrdinalType& ldz, OrdinalType* M, MagnitudeType *PL, MagnitudeType *PR, MagnitudeType *DIF, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK,
const OrdinalType& liwork, OrdinalType* info )
const;
404 void GGES(
const char& JOBVL,
const char& JOBVR,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType* , ScalarType* , ScalarType* ),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* sdim, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* BWORK, OrdinalType* info )
const;
411 void GESVD(
const char& JOBU,
const char& JOBVT,
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* S, ScalarType* U,
const OrdinalType& ldu, ScalarType* V,
const OrdinalType& ldv, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const;
428 void ORMQR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
435 void ORM2R(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k,
const ScalarType A[],
const OrdinalType& lda,
const ScalarType TAU[], ScalarType C[],
const OrdinalType& ldc, ScalarType WORK[], OrdinalType*
const info)
const;
445 void UNMQR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
451 void UNM2R(
const char& SIDE,
const char&
TRANS,
const OrdinalType& M,
const OrdinalType& N,
const OrdinalType& K,
const ScalarType A[],
const OrdinalType& LDA,
const ScalarType TAU[], ScalarType C[],
const OrdinalType& LDC, ScalarType WORK[], OrdinalType*
const INFO)
const;
462 void ORGQR(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
472 void UNGQR(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
477 void ORGHR(
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
482 void ORMHR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const;
490 void TREVC(
const char& SIDE,
const char& HOWMNY, OrdinalType* select,
const OrdinalType& n,
const ScalarType* T,
const OrdinalType& ldt, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* m, ScalarType* WORK, OrdinalType* info)
const;
495 void TREVC(
const char& SIDE,
const OrdinalType& n,
const ScalarType* T,
const OrdinalType& ldt, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* m, ScalarType* WORK, MagnitudeType* RWORK, OrdinalType* info)
const;
501 TEUCHOS_DEPRECATED
void TREXC(
const char& COMPQ,
const OrdinalType& n, ScalarType* T,
const OrdinalType& ldt, ScalarType* Q,
const OrdinalType& ldq,
const OrdinalType& ifst,
const OrdinalType& ilst, ScalarType* WORK, OrdinalType* info)
const;
502 void TREXC(
const char& COMPQ,
const OrdinalType& n, ScalarType* T,
const OrdinalType& ldt, ScalarType* Q,
const OrdinalType& ldq, OrdinalType* ifst, OrdinalType* ilst, ScalarType* WORK, OrdinalType* info)
const;
507 void TGEVC(
const char& SIDE,
const char& HOWMNY,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* S,
const OrdinalType& lds, ScalarType* P,
const OrdinalType& ldp, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* M, ScalarType* WORK, OrdinalType* info)
const;
516 void LARTG(
const ScalarType& f,
const ScalarType& g, MagnitudeType* c, ScalarType* s, ScalarType* r )
const;
519 void LARFG(
const OrdinalType& n, ScalarType* alpha, ScalarType* x,
const OrdinalType& incx, ScalarType* tau )
const;
529 TEUCHOS_DEPRECATED
void GEBAL(
const char& JOBZ,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda,
const OrdinalType& ilo,
const OrdinalType& ihi, MagnitudeType* scale, OrdinalType* info)
const;
530 void GEBAL(
const char& JOBZ,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* scale, OrdinalType* info)
const;
533 void GEBAK(
const char& JOBZ,
const char& SIDE,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi,
const MagnitudeType* scale ,
const OrdinalType& m, ScalarType* V,
const OrdinalType& ldv, OrdinalType* info)
const;
539 ScalarType
LARND(
const OrdinalType& idist, OrdinalType* seed )
const;
543 void LARNV(
const OrdinalType& idist, OrdinalType* seed,
const OrdinalType& n, ScalarType* v )
const;
551 ScalarType
LAMCH(
const char& CMACH)
const;
557 OrdinalType
ILAENV(
const OrdinalType& ispec,
const std::string& NAME,
const std::string& OPTS,
const OrdinalType& N1 = -1,
const OrdinalType& N2 = -1,
const OrdinalType& N3 = -1,
const OrdinalType& N4 = -1 )
const;
565 ScalarType
LAPY2(
const ScalarType& x,
const ScalarType& y)
const;
574 template<
typename OrdinalType,
typename ScalarType>
577 UndefinedLAPACKRoutine<ScalarType>::notDefined();
580 template<
typename OrdinalType,
typename ScalarType>
581 void LAPACK<OrdinalType, ScalarType>::PTTRS(
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* d,
const ScalarType* e, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
583 UndefinedLAPACKRoutine<ScalarType>::notDefined();
586 template<
typename OrdinalType,
typename ScalarType>
589 UndefinedLAPACKRoutine<ScalarType>::notDefined();
592 template<
typename OrdinalType,
typename ScalarType>
593 void LAPACK<OrdinalType, ScalarType>::POTRS(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
595 UndefinedLAPACKRoutine<ScalarType>::notDefined();
598 template<
typename OrdinalType,
typename ScalarType>
601 UndefinedLAPACKRoutine<ScalarType>::notDefined();
604 template<
typename OrdinalType,
typename ScalarType>
605 void LAPACK<OrdinalType, ScalarType>::POCON(
const char& UPLO,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
607 UndefinedLAPACKRoutine<ScalarType>::notDefined();
610 template<
typename OrdinalType,
typename ScalarType>
611 void LAPACK<OrdinalType, ScalarType>::POSV(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
613 UndefinedLAPACKRoutine<ScalarType>::notDefined();
616 template<
typename OrdinalType,
typename ScalarType>
617 void LAPACK<OrdinalType, ScalarType>::POEQU(
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda, MagnitudeType* S, MagnitudeType* scond, MagnitudeType* amax, OrdinalType* info)
const
623 }
else if (lda < TEUCHOS_MAX(1, n)) {
644 MagnitudeType smin = S[0];
646 for (OrdinalType i=0; i<n; ++i) {
648 smin = TEUCHOS_MIN( smin, S[i] );
649 *amax = TEUCHOS_MAX( *amax, S[i] );
654 for (OrdinalType i=0; i<n; ++i) {
660 for (OrdinalType i=0; i<n; ++i) {
668 template<
typename OrdinalType,
typename ScalarType>
669 void LAPACK<OrdinalType, ScalarType>::PORFS(
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
671 UndefinedLAPACKRoutine<ScalarType>::notDefined();
675 template<
typename OrdinalType,
typename ScalarType>
676 TEUCHOS_DEPRECATED
void LAPACK<OrdinalType, ScalarType>::POSVX(
const char& FACT,
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf,
const char& EQUED, ScalarType* S, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
678 UndefinedLAPACKRoutine<ScalarType>::notDefined();
680 template<
typename OrdinalType,
typename ScalarType>
681 void LAPACK<OrdinalType, ScalarType>::POSVX(
const char& FACT,
const char& UPLO,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf,
char* EQUED, ScalarType* S, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
683 UndefinedLAPACKRoutine<ScalarType>::notDefined();
686 template<
typename OrdinalType,
typename ScalarType>
687 void LAPACK<OrdinalType,ScalarType>::GELS(
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
689 UndefinedLAPACKRoutine<ScalarType>::notDefined();
692 template<
typename OrdinalType,
typename ScalarType>
693 void LAPACK<OrdinalType, ScalarType>::GELSS(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* S,
const MagnitudeType rcond, OrdinalType* rank, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const
695 UndefinedLAPACKRoutine<ScalarType>::notDefined();
698 template<
typename OrdinalType,
typename ScalarType>
699 void LAPACK<OrdinalType,ScalarType>::GELSS(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* S,
const ScalarType& rcond, OrdinalType* rank, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
701 UndefinedLAPACKRoutine<ScalarType>::notDefined();
704 template<
typename OrdinalType,
typename ScalarType>
705 void LAPACK<OrdinalType,ScalarType>::GGLSE(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& p, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* C, ScalarType* D, ScalarType* X, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
707 UndefinedLAPACKRoutine<ScalarType>::notDefined();
710 template<
typename OrdinalType,
typename ScalarType>
711 void LAPACK<OrdinalType,ScalarType>::GEQRF(
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
713 UndefinedLAPACKRoutine<ScalarType>::notDefined();
716 template<
typename OrdinalType,
typename ScalarType>
717 void LAPACK<OrdinalType,ScalarType>::GEQR2 (
const OrdinalType& m,
const OrdinalType& n, ScalarType A[],
const OrdinalType& lda, ScalarType TAU[], ScalarType WORK[], OrdinalType*
const info)
const
719 UndefinedLAPACKRoutine<ScalarType>::notDefined();
722 template<
typename OrdinalType,
typename ScalarType>
725 UndefinedLAPACKRoutine<ScalarType>::notDefined();
728 template<
typename OrdinalType,
typename ScalarType>
729 void LAPACK<OrdinalType,ScalarType>::GETRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
731 UndefinedLAPACKRoutine<ScalarType>::notDefined();
734 template<
typename OrdinalType,
typename ScalarType>
735 void LAPACK<OrdinalType,ScalarType>::LASCL(
const char& TYPE,
const OrdinalType& kl,
const OrdinalType& ku,
const MagnitudeType cfrom,
const MagnitudeType cto,
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* info)
const
751 MagnitudeType cfromc = cfrom;
752 MagnitudeType ctoc = cto;
753 MagnitudeType cfrom1;
758 cfrom1 = cfromc*smlnum;
759 if (cfrom1 == cfromc) {
765 cto1 = ctoc / bignum;
785 for (j=0; j<n; j++) {
787 for (i=0; i<m; i++) { *ptr = mul * (*ptr); ptr++; }
793 template<
typename OrdinalType,
typename ScalarType>
796 GEQP3 (
const OrdinalType& m,
797 const OrdinalType& n,
799 const OrdinalType& lda,
803 const OrdinalType& lwork,
804 MagnitudeType* RWORK,
805 OrdinalType* info)
const
807 UndefinedLAPACKRoutine<ScalarType>::notDefined();
810 template<
typename OrdinalType,
typename ScalarType>
813 LASWP (
const OrdinalType& N,
815 const OrdinalType& LDA,
816 const OrdinalType& K1,
817 const OrdinalType& K2,
818 const OrdinalType IPIV[],
819 const OrdinalType& INCX)
const
821 UndefinedLAPACKRoutine<ScalarType>::notDefined();
824 template<
typename OrdinalType,
typename ScalarType>
825 void LAPACK<OrdinalType,ScalarType>::GBTRF(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku, ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV, OrdinalType* info)
const
827 UndefinedLAPACKRoutine<ScalarType>::notDefined();
830 template<
typename OrdinalType,
typename ScalarType>
831 void LAPACK<OrdinalType,ScalarType>::GBTRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
833 UndefinedLAPACKRoutine<ScalarType>::notDefined();
836 template<
typename OrdinalType,
typename ScalarType>
839 UndefinedLAPACKRoutine<ScalarType>::notDefined();
842 template<
typename OrdinalType,
typename ScalarType>
843 void LAPACK<OrdinalType,ScalarType>::GTTRS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* dl,
const ScalarType* d,
const ScalarType* du,
const ScalarType* du2,
const OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
845 UndefinedLAPACKRoutine<ScalarType>::notDefined();
848 template<
typename OrdinalType,
typename ScalarType>
849 void LAPACK<OrdinalType,ScalarType>::GETRI(
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda,
const OrdinalType* IPIV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
851 UndefinedLAPACKRoutine<ScalarType>::notDefined();
854 template<
typename OrdinalType,
typename ScalarType>
857 LATRS (
const char& UPLO,
861 const OrdinalType& N,
863 const OrdinalType& LDA,
865 MagnitudeType* SCALE,
866 MagnitudeType* CNORM,
867 OrdinalType* INFO)
const
869 UndefinedLAPACKRoutine<ScalarType>::notDefined();
872 template<
typename OrdinalType,
typename ScalarType>
873 void LAPACK<OrdinalType,ScalarType>::GECON(
const char& NORM,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
875 UndefinedLAPACKRoutine<ScalarType>::notDefined();
878 template<
typename OrdinalType,
typename ScalarType>
879 void LAPACK<OrdinalType,ScalarType>::GBCON(
const char& NORM,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV,
const ScalarType& anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
881 UndefinedLAPACKRoutine<ScalarType>::notDefined();
884 template<
typename OrdinalType,
typename ScalarType>
887 UndefinedLAPACKRoutine<ScalarType>::notDefined();
890 template<
typename OrdinalType,
typename ScalarType>
891 void LAPACK<OrdinalType,ScalarType>::GESV(
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, OrdinalType* IPIV, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
893 UndefinedLAPACKRoutine<ScalarType>::notDefined();
896 template<
typename OrdinalType,
typename ScalarType>
897 void LAPACK<OrdinalType,ScalarType>::GEEQU(
const OrdinalType& m,
const OrdinalType& n,
const ScalarType* A,
const OrdinalType& lda, ScalarType* R, ScalarType* C, ScalarType* rowcond, ScalarType* colcond, ScalarType* amax, OrdinalType* info)
const
906 }
else if (lda < TEUCHOS_MAX(1, m)) {
919 if (m == 0 || n == 0) {
931 for (OrdinalType i=0; i<m; i++) {
936 for (OrdinalType j=0; j<n; j++) {
937 for (OrdinalType i=0; i<m; i++) {
943 MagnitudeType rcmin = bignum;
944 MagnitudeType rcmax = mZero;
945 for (OrdinalType i=0; i<m; i++) {
946 rcmax = TEUCHOS_MAX( rcmax, R[i] );
947 rcmin = TEUCHOS_MIN( rcmin, R[i] );
951 if (rcmin == mZero) {
953 for (OrdinalType i=0; i<m; i++) {
959 for (OrdinalType i=0; i<m; i++) {
960 R[i] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( R[i], smlnum ), bignum );
963 *rowcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum );
967 for (OrdinalType j=0; j<n; j++) {
972 for (OrdinalType j=0; j<n; j++) {
973 for (OrdinalType i=0; i<m; i++) {
981 for (OrdinalType j=0; j<n; j++) {
982 rcmax = TEUCHOS_MAX( rcmax, C[j] );
983 rcmin = TEUCHOS_MIN( rcmin, C[j] );
986 if (rcmin == mZero) {
988 for (OrdinalType j=0; j<n; j++) {
994 for (OrdinalType j=0; j<n; j++) {
995 C[j] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( C[j], smlnum ), bignum );
998 *colcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum );
1002 template<
typename OrdinalType,
typename ScalarType>
1003 void LAPACK<OrdinalType,ScalarType>::GERFS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const OrdinalType* IPIV,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
1005 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1008 template<
typename OrdinalType,
typename ScalarType>
1009 void LAPACK<OrdinalType,ScalarType>::GBEQU(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const ScalarType* A,
const OrdinalType& lda, MagnitudeType* R, MagnitudeType* C, MagnitudeType* rowcond, MagnitudeType* colcond, MagnitudeType* amax, OrdinalType* info)
const
1018 }
else if (kl < 0) {
1020 }
else if (ku < 0) {
1022 }
else if (lda < kl+ku+1) {
1035 if (m == 0 || n == 0) {
1047 for (OrdinalType i=0; i<m; i++) {
1052 for (OrdinalType j=0; j<n; j++) {
1053 for (OrdinalType i=TEUCHOS_MAX(j-ku,0); i<TEUCHOS_MIN(j+kl,m-1); i++) {
1059 MagnitudeType rcmin = bignum;
1060 MagnitudeType rcmax = mZero;
1061 for (OrdinalType i=0; i<m; i++) {
1062 rcmax = TEUCHOS_MAX( rcmax, R[i] );
1063 rcmin = TEUCHOS_MIN( rcmin, R[i] );
1067 if (rcmin == mZero) {
1069 for (OrdinalType i=0; i<m; i++) {
1075 for (OrdinalType i=0; i<m; i++) {
1076 R[i] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( R[i], smlnum ), bignum );
1079 *rowcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum );
1083 for (OrdinalType j=0; j<n; j++) {
1088 for (OrdinalType j=0; j<n; j++) {
1089 for (OrdinalType i=TEUCHOS_MAX(j-ku,0); i<TEUCHOS_MIN(j+kl,m-1); i++) {
1097 for (OrdinalType j=0; j<n; j++) {
1098 rcmax = TEUCHOS_MAX( rcmax, C[j] );
1099 rcmin = TEUCHOS_MIN( rcmin, C[j] );
1102 if (rcmin == mZero) {
1104 for (OrdinalType j=0; j<n; j++) {
1110 for (OrdinalType j=0; j<n; j++) {
1111 C[j] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( C[j], smlnum ), bignum );
1114 *colcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum );
1118 template<
typename OrdinalType,
typename ScalarType>
1119 void LAPACK<OrdinalType,ScalarType>::GBRFS(
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& kl,
const OrdinalType& ku,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* AF,
const OrdinalType& ldaf,
const OrdinalType* IPIV,
const ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
1121 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1125 template<
typename OrdinalType,
typename ScalarType>
1126 TEUCHOS_DEPRECATED
void LAPACK<OrdinalType,ScalarType>::GESVX(
const char& FACT,
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf, OrdinalType* IPIV,
const char& EQUED, ScalarType* R, ScalarType* C, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
1128 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1130 template<
typename OrdinalType,
typename ScalarType>
1131 void LAPACK<OrdinalType,ScalarType>::GESVX(
const char& FACT,
const char&
TRANS,
const OrdinalType& n,
const OrdinalType& nrhs, ScalarType* A,
const OrdinalType& lda, ScalarType* AF,
const OrdinalType& ldaf, OrdinalType* IPIV,
char* EQUED, ScalarType* R, ScalarType* C, ScalarType* B,
const OrdinalType& ldb, ScalarType* X,
const OrdinalType& ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
1133 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1136 template<
typename OrdinalType,
typename ScalarType>
1137 void LAPACK<OrdinalType,ScalarType>::SYTRD(
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* D, ScalarType* E, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1139 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1142 template<
typename OrdinalType,
typename ScalarType>
1143 void LAPACK<OrdinalType,ScalarType>::GEHRD(
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* A,
const OrdinalType& lda, ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1145 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1148 template<
typename OrdinalType,
typename ScalarType>
1149 void LAPACK<OrdinalType,ScalarType>::TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const OrdinalType& n,
const OrdinalType& nrhs,
const ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* info)
const
1151 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1154 template<
typename OrdinalType,
typename ScalarType>
1157 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1160 template<
typename OrdinalType,
typename ScalarType>
1161 void LAPACK<OrdinalType,ScalarType>::SPEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* AP, ScalarType* W, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK, OrdinalType* info)
const
1163 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1166 template<
typename OrdinalType,
typename ScalarType>
1167 void LAPACK<OrdinalType,ScalarType>::SYEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* W, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1169 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1172 template<
typename OrdinalType,
typename ScalarType>
1173 void LAPACK<OrdinalType,ScalarType>::SYGV(
const OrdinalType& itype,
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, ScalarType* W, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1175 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1178 template<
typename OrdinalType,
typename ScalarType>
1179 void LAPACK<OrdinalType,ScalarType>::HEEV(
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* W, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const
1181 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1184 template<
typename OrdinalType,
typename ScalarType>
1185 void LAPACK<OrdinalType,ScalarType>::HEGV(
const OrdinalType& itype,
const char& JOBZ,
const char& UPLO,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* W, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const
1187 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1190 template<
typename OrdinalType,
typename ScalarType>
1191 void LAPACK<OrdinalType,ScalarType>::STEQR(
const char& COMPZ,
const OrdinalType& n, ScalarType* D, ScalarType* E, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK, OrdinalType* info)
const
1193 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1196 template<
typename OrdinalType,
typename ScalarType>
1197 void LAPACK<OrdinalType, ScalarType>::HSEQR(
const char& JOB,
const char& COMPZ,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* H,
const OrdinalType& ldh, ScalarType* WR, ScalarType* WI, ScalarType* Z,
const OrdinalType& ldz, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1199 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1202 template<
typename OrdinalType,
typename ScalarType>
1203 void LAPACK<OrdinalType, ScalarType>::GEES(
const char& JOBVS,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType*, ScalarType*),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, ScalarType* WR, ScalarType* WI, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* BWORK, OrdinalType* info)
const
1205 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1208 template<
typename OrdinalType,
typename ScalarType>
1209 void LAPACK<OrdinalType, ScalarType>::GEES(
const char& JOBVS,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType*),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, ScalarType* W, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType *RWORK, OrdinalType* BWORK, OrdinalType* info)
const
1211 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1214 template<
typename OrdinalType,
typename ScalarType>
1215 void LAPACK<OrdinalType, ScalarType>::GEES(
const char& JOBVS,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* sdim, MagnitudeType* WR, MagnitudeType* WI, ScalarType* VS,
const OrdinalType& ldvs, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType *RWORK, OrdinalType* BWORK, OrdinalType* info)
const
1217 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1220 template<
typename OrdinalType,
typename ScalarType>
1221 void LAPACK<OrdinalType, ScalarType>::GEEV(
const char& JOBVL,
const char& JOBVR,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* WR, MagnitudeType* WI, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* rwork, OrdinalType* info)
const
1223 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1226 template<
typename OrdinalType,
typename ScalarType>
1227 void LAPACK<OrdinalType, ScalarType>::GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* WR, ScalarType* WI, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* SCALE, MagnitudeType* abnrm, MagnitudeType* RCONDE, MagnitudeType* RCONDV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK, OrdinalType* info)
const
1229 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1232 template<
typename OrdinalType,
typename ScalarType>
1233 void LAPACK<OrdinalType, ScalarType>::GESVD(
const char& JOBU,
const char& JOBVT,
const OrdinalType& m,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, MagnitudeType* S, ScalarType* U,
const OrdinalType& ldu, ScalarType* V,
const OrdinalType& ldv, ScalarType* WORK,
const OrdinalType& lwork, MagnitudeType* RWORK, OrdinalType* info)
const
1235 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1238 template<
typename OrdinalType,
typename ScalarType>
1239 void LAPACK<OrdinalType, ScalarType>::GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType* ALPHAR, MagnitudeType* ALPHAI, ScalarType* BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* lscale, MagnitudeType* rscale, MagnitudeType* abnrm, MagnitudeType* bbnrm, MagnitudeType* RCONDE, MagnitudeType* RCONDV, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK, OrdinalType* BWORK, OrdinalType* info)
const
1241 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1244 template<
typename OrdinalType,
typename ScalarType>
1245 void LAPACK<OrdinalType, ScalarType>::GGEV(
const char& JOBVL,
const char& JOBVR,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType* BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1247 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1251 template<
typename OrdinalType,
typename ScalarType>
1252 void LAPACK<OrdinalType,ScalarType>::TRSEN(
const char& JOB,
const char& COMPQ,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* T,
const OrdinalType& ldt, ScalarType* Q,
const OrdinalType& ldq, MagnitudeType *WR, MagnitudeType *WI, OrdinalType* M, ScalarType* S, MagnitudeType *SEP, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK,
const OrdinalType& liwork, OrdinalType* info )
const
1254 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1258 template<
typename OrdinalType,
typename ScalarType>
1259 void LAPACK<OrdinalType,ScalarType>::TGSEN(
const OrdinalType& ijob,
const OrdinalType& wantq,
const OrdinalType& wantz,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType* Q,
const OrdinalType& ldq, ScalarType* Z,
const OrdinalType& ldz, OrdinalType* M, MagnitudeType *PL, MagnitudeType *PR, MagnitudeType *DIF, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* IWORK,
const OrdinalType& liwork, OrdinalType* info )
const
1261 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1265 template<
typename OrdinalType,
typename ScalarType>
1266 void LAPACK<OrdinalType, ScalarType>::GGES(
const char& JOBVL,
const char& JOBVR,
const char& SORT, OrdinalType& (*ptr2func)(ScalarType*, ScalarType*, ScalarType*),
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, ScalarType* B,
const OrdinalType& ldb, OrdinalType* sdim, MagnitudeType* ALPHAR, MagnitudeType* ALPHAI, MagnitudeType* BETA, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* BWORK, OrdinalType* info )
const
1268 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1271 template<
typename OrdinalType,
typename ScalarType>
1272 void LAPACK<OrdinalType, ScalarType>::ORMQR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1274 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1277 template<
typename OrdinalType,
typename ScalarType>
1278 void LAPACK<OrdinalType, ScalarType>::ORM2R(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k,
const ScalarType A[],
const OrdinalType& lda,
const ScalarType TAU[], ScalarType C[],
const OrdinalType& ldc, ScalarType WORK[], OrdinalType*
const info)
const
1280 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1283 template<
typename OrdinalType,
typename ScalarType>
1284 void LAPACK<OrdinalType, ScalarType>::UNMQR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1286 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1289 template<
typename OrdinalType,
typename ScalarType>
1290 void LAPACK<OrdinalType, ScalarType>::UNM2R(
const char& SIDE,
const char&
TRANS,
const OrdinalType& M,
const OrdinalType& N,
const OrdinalType& K,
const ScalarType A[],
const OrdinalType& LDA,
const ScalarType TAU[], ScalarType C[],
const OrdinalType& LDC, ScalarType WORK[], OrdinalType*
const INFO)
const
1292 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1295 template<
typename OrdinalType,
typename ScalarType>
1296 void LAPACK<OrdinalType, ScalarType>::ORGQR(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1298 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1301 template<
typename OrdinalType,
typename ScalarType>
1302 void LAPACK<OrdinalType, ScalarType>::UNGQR(
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& k, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1304 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1307 template<
typename OrdinalType,
typename ScalarType>
1308 void LAPACK<OrdinalType, ScalarType>::ORGHR(
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi, ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1310 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1313 template<
typename OrdinalType,
typename ScalarType>
1314 void LAPACK<OrdinalType, ScalarType>::ORMHR(
const char& SIDE,
const char&
TRANS,
const OrdinalType& m,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi,
const ScalarType* A,
const OrdinalType& lda,
const ScalarType* TAU, ScalarType* C,
const OrdinalType& ldc, ScalarType* WORK,
const OrdinalType& lwork, OrdinalType* info)
const
1316 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1319 template<
typename OrdinalType,
typename ScalarType>
1320 void LAPACK<OrdinalType, ScalarType>::TREVC(
const char& SIDE,
const char& HOWMNY, OrdinalType* select,
const OrdinalType& n,
const ScalarType* T,
const OrdinalType& ldt, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* m, ScalarType* WORK, OrdinalType* info)
const
1322 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1325 template<
typename OrdinalType,
typename ScalarType>
1326 void LAPACK<OrdinalType, ScalarType>::TREVC(
const char& SIDE,
const OrdinalType& n,
const ScalarType* T,
const OrdinalType& ldt, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* m, ScalarType* WORK, MagnitudeType* RWORK, OrdinalType* info)
const
1328 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1331 template<
typename OrdinalType,
typename ScalarType>
1332 void LAPACK<OrdinalType, ScalarType>::TREXC(
const char& COMPQ,
const OrdinalType& n, ScalarType* T,
const OrdinalType& ldt, ScalarType* Q,
const OrdinalType& ldq, OrdinalType* ifst, OrdinalType* ilst, ScalarType* WORK, OrdinalType* info)
const
1334 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1338 template<
typename OrdinalType,
typename ScalarType>
1339 void LAPACK<OrdinalType, ScalarType>::TGEVC(
const char& SIDE,
const char& HOWMNY,
const OrdinalType* SELECT,
const OrdinalType& n, ScalarType* S,
const OrdinalType& lds, ScalarType* P,
const OrdinalType& ldp, ScalarType* VL,
const OrdinalType& ldvl, ScalarType* VR,
const OrdinalType& ldvr,
const OrdinalType& mm, OrdinalType* M, ScalarType* WORK, OrdinalType* info)
const
1341 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1345 template<
typename OrdinalType,
typename ScalarType>
1348 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
1351 template<
typename OrdinalType,
typename ScalarType>
1352 OrdinalType
LAPACK<OrdinalType, ScalarType>::ILAENV(
const OrdinalType& ispec,
const std::string& NAME,
const std::string& OPTS,
const OrdinalType& N1,
const OrdinalType& N2,
const OrdinalType& N3,
const OrdinalType& N4 )
const
1354 return UndefinedLAPACKRoutine<OrdinalType>::notDefined();
1357 template<
typename OrdinalType,
typename ScalarType>
1360 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
1363 template<
typename OrdinalType,
typename ScalarType>
1366 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1369 template<
typename OrdinalType,
typename ScalarType>
1372 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1375 template<
typename OrdinalType,
typename ScalarType>
1376 void LAPACK<OrdinalType, ScalarType>::GEBAL(
const char& JOBZ,
const OrdinalType& n, ScalarType* A,
const OrdinalType& lda, OrdinalType* ilo, OrdinalType* ihi, MagnitudeType* scale, OrdinalType* info )
const
1378 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1382 template<
typename OrdinalType,
typename ScalarType>
1383 void LAPACK<OrdinalType, ScalarType>::GEBAK(
const char& JOBZ,
const char& SIDE,
const OrdinalType& n,
const OrdinalType& ilo,
const OrdinalType& ihi,
const MagnitudeType* scale,
const OrdinalType& m, ScalarType* V,
const OrdinalType& ldv, OrdinalType* info )
const
1385 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1388 template<
typename OrdinalType,
typename ScalarType>
1391 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
1394 template<
typename OrdinalType,
typename ScalarType>
1397 UndefinedLAPACKRoutine<ScalarType>::notDefined();
1402 #ifndef DOXYGEN_SHOULD_SKIP_THIS
1407 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, float>
1410 inline LAPACK(
void) {}
1411 inline LAPACK(
const LAPACK<int, float>& ) {}
1412 inline virtual ~LAPACK(
void) {}
1415 void POTRF(
const char& UPLO,
const int& n,
float* A,
const int& lda,
int* info)
const;
1416 void POTRS(
const char& UPLO,
const int& n,
const int& nrhs,
const float* A,
const int& lda,
float* B,
const int& ldb,
int* info)
const;
1417 void PTTRF(
const int& n,
float* d,
float* e,
int* info)
const;
1418 void PTTRS(
const int& n,
const int& nrhs,
const float* d,
const float* e,
float* B,
const int& ldb,
int* info)
const;
1419 void POTRI(
const char& UPLO,
const int& n,
float* A,
const int& lda,
int* info)
const;
1420 void POCON(
const char& UPLO,
const int& n,
const float* A,
const int& lda,
const float& anorm,
float* rcond,
float* WORK,
int* IWORK,
int* info)
const;
1421 void POSV(
const char& UPLO,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* B,
const int& ldb,
int* info)
const;
1422 void POEQU(
const int& n,
const float* A,
const int& lda,
float* S,
float* scond,
float* amax,
int* info)
const;
1423 void PORFS(
const char& UPLO,
const int& n,
const int& nrhs,
float* A,
const int& lda,
const float* AF,
const int& ldaf,
const float* B,
const int& ldb,
float* X,
const int& ldx,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1426 TEUCHOS_DEPRECATED
void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* AF,
const int& ldaf,
const char& EQUED,
float* S,
float* B,
const int& ldb,
float* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1427 void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* AF,
const int& ldaf,
char* EQUED,
float* S,
float* B,
const int& ldb,
float* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1430 void GELS(
const char&
TRANS,
const int& m,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* B,
const int& ldb,
float* WORK,
const int& lwork,
int* info)
const;
1431 void GELSS(
const int& m,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* B,
const int& ldb,
float* S,
const float& rcond,
int* rank,
float* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1432 void GELSS(
const int& m,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* B,
const int& ldb,
float* S,
const float& rcond,
int* rank,
float* WORK,
const int& lwork,
int* info)
const;
1433 void GGLSE(
const int& m,
const int& n,
const int& p,
float* A,
const int& lda,
float* B,
const int& ldb,
float* C,
float* D,
float* X,
float* WORK,
const int& lwork,
int* info)
const;
1434 void GEQRF(
const int& m,
const int& n,
float* A,
const int& lda,
float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1435 void GEQR2(
const int& m,
const int& n,
float A[],
const int& lda,
float TAU[],
float WORK[],
int*
const info)
const;
1437 void GETRF(
const int& m,
const int& n,
float* A,
const int& lda,
int* IPIV,
int* info)
const;
1438 void GETRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const float* A,
const int& lda,
const int* IPIV,
float* B,
const int& ldb,
int* info)
const;
1439 void LASCL(
const char& TYPE,
const int& kl,
const int& ku,
const float& cfrom,
const float& cto,
const int& m,
const int& n,
float* A,
const int& lda,
int* info)
const;
1442 GEQP3 (
const int& m,
1453 void LASWP (
const int& N,
1459 const int& INCX)
const;
1461 void GBTRF(
const int& m,
const int& n,
const int& kl,
const int& ku,
float* A,
const int& lda,
int* IPIV,
int* info)
const;
1462 void GBTRS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const float* A,
const int& lda,
const int* IPIV,
float* B,
const int& ldb,
int* info)
const;
1463 void GTTRF(
const int& n,
float* dl,
float* d,
float* du,
float* du2,
int* IPIV,
int* info)
const;
1464 void GTTRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const float* dl,
const float* d,
const float* du,
const float* du2,
const int* IPIV,
float* B,
const int& ldb,
int* info)
const;
1467 void GETRI(
const int& n,
float* A,
const int& lda,
const int* IPIV,
float* WORK,
const int& lwork,
int* info)
const;
1468 void LATRS (
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const char& NORMIN,
const int& N,
float* A,
const int& LDA,
float* X,
float* SCALE,
float* CNORM,
int* INFO)
const;
1469 void GECON(
const char& NORM,
const int& n,
const float* A,
const int& lda,
const float& anorm,
float* rcond,
float* WORK,
int* IWORK,
int* info)
const;
1470 void GBCON(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const float* A,
const int& lda,
int* IPIV,
const float& anorm,
float* rcond,
float* WORK,
int* IWORK,
int* info)
const;
1471 float LANGB(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const float* A,
const int& lda,
float* WORK)
const;
1472 void GESV(
const int& n,
const int& nrhs,
float* A,
const int& lda,
int* IPIV,
float* B,
const int& ldb,
int* info)
const;
1473 void GEEQU(
const int& m,
const int& n,
const float* A,
const int& lda,
float* R,
float* C,
float* rowcond,
float* colcond,
float* amax,
int* info)
const;
1474 void GERFS(
const char&
TRANS,
const int& n,
const int& nrhs,
const float* A,
const int& lda,
const float* AF,
const int& ldaf,
const int* IPIV,
const float* B,
const int& ldb,
float* X,
const int& ldx,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1475 void GBEQU(
const int& m,
const int& n,
const int& kl,
const int& ku,
const float* A,
const int& lda,
float* R,
float* C,
float* rowcond,
float* colcond,
float* amax,
int* info)
const;
1476 void GBRFS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const float* A,
const int& lda,
const float* AF,
const int& ldaf,
const int* IPIV,
const float* B,
const int& ldb,
float* X,
const int& ldx,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1479 TEUCHOS_DEPRECATED
void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* AF,
const int& ldaf,
int* IPIV,
const char& EQUED,
float* R,
float* C,
float* B,
const int& ldb,
float* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1480 void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs,
float* A,
const int& lda,
float* AF,
const int& ldaf,
int* IPIV,
char* EQUED,
float* R,
float* C,
float* B,
const int& ldb,
float* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR,
float* WORK,
int* IWORK,
int* info)
const;
1482 void SYTRD(
const char& UPLO,
const int& n,
float* A,
const int& lda,
float* D,
float* E,
float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1483 void GEHRD(
const int& n,
const int& ilo,
const int& ihi,
float* A,
const int& lda,
float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1484 void TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const int& n,
const int& nrhs,
const float* A,
const int& lda,
float* B,
const int& ldb,
int* info)
const;
1485 void TRTRI(
const char& UPLO,
const char& DIAG,
const int& n,
const float* A,
const int& lda,
int* info)
const;
1488 void SPEV(
const char& JOBZ,
const char& UPLO,
const int& n,
float* AP,
float* W,
float* Z,
const int& ldz,
float* WORK,
int* info)
const;
1489 void SYEV(
const char& JOBZ,
const char& UPLO,
const int& n,
float* A,
const int& lda,
float* W,
float* WORK,
const int& lwork,
int* info)
const;
1490 void SYGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* W,
float* WORK,
const int& lwork,
int* info)
const;
1491 void HEEV(
const char& JOBZ,
const char& UPLO,
const int& n,
float* A,
const int& lda,
float* W,
float* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1492 void HEGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* W,
float* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1493 void STEQR(
const char& COMPZ,
const int& n,
float* D,
float* E,
float* Z,
const int& ldz,
float* WORK,
int* info)
const;
1496 void HSEQR(
const char& JOB,
const char& COMPZ,
const int& n,
const int& ilo,
const int& ihi,
float* H,
const int& ldh,
float* WR,
float* WI,
float* Z,
const int& ldz,
float* WORK,
const int& lwork,
int* info)
const;
1497 void GEES(
const char& JOBVS,
const char& SORT,
int (*ptr2func)(
float*,
float*),
const int& n,
float* A,
const int& lda,
int* sdim,
float* WR,
float* WI,
float* VS,
const int& ldvs,
float* WORK,
const int& lwork,
int* BWORK,
int* info)
const;
1498 void GEES(
const char& JOBVS,
const int& n,
float* A,
const int& lda,
int* sdim,
float* WR,
float* WI,
float* VS,
const int& ldvs,
float* WORK,
const int& lwork,
float* RWORK,
int* BWORK,
int* info)
const;
1500 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
float* A,
const int& lda,
float* WR,
float* WI,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
float* WORK,
const int& lwork,
int* info)
const;
1501 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
float* A,
const int& lda,
float* WR,
float* WI,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
float* WORK,
const int& lwork,
float* rwork,
int* info)
const;
1503 void GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
float* A,
const int& lda,
float* WR,
float* WI,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* SCALE,
float* abnrm,
float* RCONDE,
float* RCONDV,
float* WORK,
const int& lwork,
int* IWORK,
int* info)
const;
1504 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* ALPHAR,
float* ALPHAI,
float* BETA,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* lscale,
float* rscale,
float* abnrm,
float* bbnrm,
float* RCONDE,
float* RCONDV,
float* WORK,
const int& lwork,
int* IWORK,
int* BWORK,
int* info)
const;
1505 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* ALPHAR,
float* ALPHAI,
float* BETA,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* lscale,
float* rscale,
float* abnrm,
float* bbnrm,
float* RCONDE,
float* RCONDV,
float* WORK,
const int& lwork,
float* rwork,
int* IWORK,
int* BWORK,
int* info)
const;
1506 void GGEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* ALPHAR,
float* ALPHAI,
float* BETA,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
float* WORK,
const int& lwork,
int* info)
const;
1507 void TRSEN(
const char& JOB,
const char& COMPQ,
const int* SELECT,
const int& n,
float* T,
const int& ldt,
float* Q,
const int& ldq,
float* WR,
float* WI,
int* M,
float* S,
float* SEP,
float* WORK,
const int& lwork,
int* IWORK,
const int& liwork,
int* info )
const;
1508 void TGSEN(
const int& ijob,
const int& wantq,
const int& wantz,
const int* SELECT,
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
float* ALPHAR,
float* ALPHAI,
float* BETA,
float* Q,
const int& ldq,
float* Z,
const int& ldz,
int* M,
float* PL,
float* PR,
float* DIF,
float* WORK,
const int& lwork,
int* IWORK,
const int& liwork,
int* info )
const;
1509 void GGES(
const char& JOBVL,
const char& JOBVR,
const char& SORT,
int (*ptr2func)(
float*,
float*,
float*),
const int& n,
float* A,
const int& lda,
float* B,
const int& ldb,
int* sdim,
float* ALPHAR,
float* ALPHAI,
float* BETA,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
float* WORK,
const int& lwork,
int* bwork,
int* info )
const;
1512 void GESVD(
const char& JOBU,
const char& JOBVT,
const int& m,
const int& n,
float* A,
const int& lda,
float* S,
float* U,
const int& ldu,
float* V,
const int& ldv,
float* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1515 void ORMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
float* A,
const int& lda,
const float* TAU,
float* C,
const int& ldc,
float* WORK,
const int& lwork,
int* info)
const;
1516 void ORM2R(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
const float A[],
const int& lda,
const float TAU[],
float C[],
const int& ldc,
float WORK[],
int*
const info)
const;
1517 void UNMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
float* A,
const int& lda,
const float* TAU,
float* C,
const int& ldc,
float* WORK,
const int& lwork,
int* info)
const;
1518 void UNM2R(
const char& SIDE,
const char&
TRANS,
const int& M,
const int& N,
const int& K,
const float A[],
const int& LDA,
const float TAU[],
float C[],
const int& LDC,
float WORK[],
int*
const INFO)
const;
1519 void ORGQR(
const int& m,
const int& n,
const int& k,
float* A,
const int& lda,
const float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1520 void UNGQR(
const int& m,
const int& n,
const int& k,
float* A,
const int& lda,
const float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1521 void ORGHR(
const int& n,
const int& ilo,
const int& ihi,
float* A,
const int& lda,
const float* TAU,
float* WORK,
const int& lwork,
int* info)
const;
1522 void ORMHR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& ilo,
const int& ihi,
const float* A,
const int& lda,
const float* TAU,
float* C,
const int& ldc,
float* WORK,
const int& lwork,
int* info)
const;
1525 void TREVC(
const char& SIDE,
const char& HOWMNY,
int* select,
const int& n,
const float* T,
const int& ldt,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
const int& mm,
int* m,
float* WORK,
int* info)
const;
1526 void TREVC(
const char& SIDE,
const int& n,
const float* T,
const int& ldt,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
const int& mm,
int* m,
float* WORK,
float* RWORK,
int* info)
const;
1529 TEUCHOS_DEPRECATED
void TREXC(
const char& COMPQ,
const int& n,
float* T,
const int& ldt,
float* Q,
const int& ldq,
const int& ifst,
const int& ilst,
float* WORK,
int* info)
const;
1530 void TREXC(
const char& COMPQ,
const int& n,
float* T,
const int& ldt,
float* Q,
const int& ldq,
int* ifst,
int* ilst,
float* WORK,
int* info)
const;
1532 void TGEVC(
const char& SIDE,
const char& HOWMNY,
const int* SELECT,
const int& n,
float* S,
const int& lds,
float* P,
const int& ldp,
float* VL,
const int& ldvl,
float* VR,
const int& ldvr,
const int& mm,
int* M,
float* WORK,
int* info)
const;
1535 void LARTG(
const float& f,
const float& g,
float* c,
float* s,
float* r )
const;
1536 void LARFG(
const int& n,
float* alpha,
float* x,
const int& incx,
float* tau )
const;
1541 TEUCHOS_DEPRECATED
void GEBAL(
const char& JOBZ,
const int& n,
float* A,
const int& lda,
const int& ilo,
const int& ihi,
float* scale,
int* info)
const;
1542 void GEBAL(
const char& JOBZ,
const int& n,
float* A,
const int& lda,
int* ilo,
int* ihi,
float* scale,
int* info)
const;
1544 void GEBAK(
const char& JOBZ,
const char& SIDE,
const int& n,
const int& ilo,
const int& ihi,
const float* scale,
const int& m,
float* V,
const int& ldv,
int* info)
const;
1547 float LARND(
const int& idist,
int* seed )
const;
1548 void LARNV(
const int& idist,
int* seed,
const int& n,
float* v )
const;
1551 float LAMCH(
const char& CMACH)
const;
1552 int ILAENV(
const int& ispec,
const std::string& NAME,
const std::string& OPTS,
const int& N1 = -1,
const int& N2 = -1,
const int& N3 = -1,
const int& N4 = -1 )
const;
1555 float LAPY2(
const float& x,
const float& y)
const;
1564 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, double>
1567 inline LAPACK(
void) {}
1568 inline LAPACK(
const LAPACK<int, double>& ) {}
1569 inline virtual ~LAPACK(
void) {}
1572 void PTTRF(
const int& n,
double* d,
double* e,
int* info)
const;
1573 void PTTRS(
const int& n,
const int& nrhs,
const double* d,
const double* e,
double* B,
const int& ldb,
int* info)
const;
1574 void POTRF(
const char& UPLO,
const int& n,
double* A,
const int& lda,
int* info)
const;
1575 void POTRS(
const char& UPLO,
const int& n,
const int& nrhs,
const double* A,
const int& lda,
double* B,
const int& ldb,
int* info)
const;
1576 void POTRI(
const char& UPLO,
const int& n,
double* A,
const int& lda,
int* info)
const;
1577 void POCON(
const char& UPLO,
const int& n,
const double* A,
const int& lda,
const double& anorm,
double* rcond,
double* WORK,
int* IWORK,
int* info)
const;
1578 void POSV(
const char& UPLO,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* B,
const int& ldb,
int* info)
const;
1579 void POEQU(
const int& n,
const double* A,
const int& lda,
double* S,
double* scond,
double* amax,
int* info)
const;
1580 void PORFS(
const char& UPLO,
const int& n,
const int& nrhs,
double* A,
const int& lda,
const double* AF,
const int& ldaf,
const double* B,
const int& ldb,
double* X,
const int& ldx,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1583 TEUCHOS_DEPRECATED
void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* AF,
const int& ldaf,
const char& EQUED,
double* S,
double* B,
const int& ldb,
double* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1584 void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* AF,
const int& ldaf,
char* EQUED,
double* S,
double* B,
const int& ldb,
double* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1587 void GELS(
const char&
TRANS,
const int& m,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* B,
const int& ldb,
double* WORK,
const int& lwork,
int* info)
const;
1588 void GELSS(
const int& m,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* B,
const int& ldb,
double* S,
const double& rcond,
int* rank,
double* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1589 void GELSS(
const int& m,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* B,
const int& ldb,
double* S,
const double& rcond,
int* rank,
double* WORK,
const int& lwork,
int* info)
const;
1590 void GGLSE(
const int& m,
const int& n,
const int& p,
double* A,
const int& lda,
double* B,
const int& ldb,
double* C,
double* D,
double* X,
double* WORK,
const int& lwork,
int* info)
const;
1591 void GEQRF(
const int& m,
const int& n,
double* A,
const int& lda,
double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1592 void GEQR2(
const int& m,
const int& n,
double A[],
const int& lda,
double TAU[],
double WORK[],
int*
const info)
const;
1593 void GETRF(
const int& m,
const int& n,
double* A,
const int& lda,
int* IPIV,
int* info)
const;
1594 void GETRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const double* A,
const int& lda,
const int* IPIV,
double* B,
const int& ldb,
int* info)
const;
1595 void LASCL(
const char& TYPE,
const int& kl,
const int& ku,
const double& cfrom,
const double& cto,
const int& m,
const int& n,
double* A,
const int& lda,
int* info)
const;
1598 GEQP3 (
const int& m,
1609 void LASWP (
const int& N,
1615 const int& INCX)
const;
1617 void GBTRF(
const int& m,
const int& n,
const int& kl,
const int& ku,
double* A,
const int& lda,
int* IPIV,
int* info)
const;
1618 void GBTRS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const double* A,
const int& lda,
const int* IPIV,
double* B,
const int& ldb,
int* info)
const;
1619 void GTTRF(
const int& n,
double* dl,
double* d,
double* du,
double* du2,
int* IPIV,
int* info)
const;
1620 void GTTRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const double* dl,
const double* d,
const double* du,
const double* du2,
const int* IPIV,
double* B,
const int& ldb,
int* info)
const;
1621 void GETRI(
const int& n,
double* A,
const int& lda,
const int* IPIV,
double* WORK,
const int& lwork,
int* info)
const;
1622 void LATRS (
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const char& NORMIN,
const int& N,
double* A,
const int& LDA,
double* X,
double* SCALE,
double* CNORM,
int* INFO)
const;
1623 void GECON(
const char& NORM,
const int& n,
const double* A,
const int& lda,
const double& anorm,
double* rcond,
double* WORK,
int* IWORK,
int* info)
const;
1624 void GBCON(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const double* A,
const int& lda,
int* IPIV,
const double& anorm,
double* rcond,
double* WORK,
int* IWORK,
int* info)
const;
1625 double LANGB(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const double* A,
const int& lda,
double* WORK)
const;
1626 void GESV(
const int& n,
const int& nrhs,
double* A,
const int& lda,
int* IPIV,
double* B,
const int& ldb,
int* info)
const;
1627 void GEEQU(
const int& m,
const int& n,
const double* A,
const int& lda,
double* R,
double* C,
double* rowcond,
double* colcond,
double* amax,
int* info)
const;
1628 void GERFS(
const char&
TRANS,
const int& n,
const int& nrhs,
const double* A,
const int& lda,
const double* AF,
const int& ldaf,
const int* IPIV,
const double* B,
const int& ldb,
double* X,
const int& ldx,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1629 void GBEQU(
const int& m,
const int& n,
const int& kl,
const int& ku,
const double* A,
const int& lda,
double* R,
double* C,
double* rowcond,
double* colcond,
double* amax,
int* info)
const;
1630 void GBRFS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const double* A,
const int& lda,
const double* AF,
const int& ldaf,
const int* IPIV,
const double* B,
const int& ldb,
double* X,
const int& ldx,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1633 TEUCHOS_DEPRECATED
void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* AF,
const int& ldaf,
int* IPIV,
const char& EQUED,
double* R,
double* C,
double* B,
const int& ldb,
double* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1634 void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs,
double* A,
const int& lda,
double* AF,
const int& ldaf,
int* IPIV,
char* EQUED,
double* R,
double* C,
double* B,
const int& ldb,
double* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR,
double* WORK,
int* IWORK,
int* info)
const;
1636 void SYTRD(
const char& UPLO,
const int& n,
double* A,
const int& lda,
double* D,
double* E,
double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1637 void GEHRD(
const int& n,
const int& ilo,
const int& ihi,
double* A,
const int& lda,
double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1638 void TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const int& n,
const int& nrhs,
const double* A,
const int& lda,
double* B,
const int& ldb,
int* info)
const;
1639 void TRTRI(
const char& UPLO,
const char& DIAG,
const int& n,
const double* A,
const int& lda,
int* info)
const;
1642 void SPEV(
const char& JOBZ,
const char& UPLO,
const int& n,
double* AP,
double* W,
double* Z,
const int& ldz,
double* WORK,
int* info)
const;
1643 void SYEV(
const char& JOBZ,
const char& UPLO,
const int& n,
double* A,
const int& lda,
double* W,
double* WORK,
const int& lwork,
int* info)
const;
1644 void SYGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* W,
double* WORK,
const int& lwork,
int* info)
const;
1645 void HEEV(
const char& JOBZ,
const char& UPLO,
const int& n,
double* A,
const int& lda,
double* W,
double* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1646 void HEGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* W,
double* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1647 void STEQR(
const char& COMPZ,
const int& n,
double* D,
double* E,
double* Z,
const int& ldz,
double* WORK,
int* info)
const;
1650 void HSEQR(
const char& JOB,
const char& COMPZ,
const int& n,
const int& ilo,
const int& ihi,
double* H,
const int& ldh,
double* WR,
double* WI,
double* Z,
const int& ldz,
double* WORK,
const int& lwork,
int* info)
const;
1651 void GEES(
const char& JOBVS,
const char& SORT,
int (*ptr2func)(
double*,
double*),
const int& n,
double* A,
const int& lda,
int* sdim,
double* WR,
double* WI,
double* VS,
const int& ldvs,
double* WORK,
const int& lwork,
int* BWORK,
int* info)
const;
1652 void GEES(
const char& JOBVS,
const int& n,
double* A,
const int& lda,
int* sdim,
double* WR,
double* WI,
double* VS,
const int& ldvs,
double* WORK,
const int& lwork,
double* RWORK,
int* BWORK,
int* info)
const;
1654 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
double* A,
const int& lda,
double* WR,
double* WI,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
double* WORK,
const int& lwork,
int* info)
const;
1655 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
double* A,
const int& lda,
double* WR,
double* WI,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
double* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1657 void GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
double* A,
const int& lda,
double* WR,
double* WI,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* SCALE,
double* abnrm,
double* RCONDE,
double* RCONDV,
double* WORK,
const int& lwork,
int* IWORK,
int* info)
const;
1658 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* ALPHAR,
double* ALPHAI,
double* BETA,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* lscale,
double* rscale,
double* abnrm,
double* bbnrm,
double* RCONDE,
double* RCONDV,
double* WORK,
const int& lwork,
int* IWORK,
int* BWORK,
int* info)
const;
1659 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* ALPHAR,
double* ALPHAI,
double* BETA,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* lscale,
double* rscale,
double* abnrm,
double* bbnrm,
double* RCONDE,
double* RCONDV,
double* WORK,
const int& lwork,
double* rwork,
int* IWORK,
int* BWORK,
int* info)
const;
1660 void GGEV(
const char& JOBVL,
const char& JOBVR,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* ALPHAR,
double* ALPHAI,
double* BETA,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
double* WORK,
const int& lwork,
int* info)
const;
1661 void TRSEN(
const char& JOB,
const char& COMPQ,
const int* SELECT,
const int& n,
double* T,
const int& ldt,
double* Q,
const int& ldq,
double* WR,
double* WI,
int* M,
double* S,
double* SEP,
double* WORK,
const int& lwork,
int* IWORK,
const int& liwork,
int* info )
const;
1662 void TGSEN(
const int& ijob,
const int& wantq,
const int& wantz,
const int* SELECT,
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
double* ALPHAR,
double* ALPHAI,
double* BETA,
double* Q,
const int& ldq,
double* Z,
const int& ldz,
int* M,
double* PL,
double* PR,
double* DIF,
double* WORK,
const int& lwork,
int* IWORK,
const int& liwork,
int* info )
const;
1663 void GGES(
const char& JOBVL,
const char& JOBVR,
const char& SORT,
int (*ptr2func)(
double*,
double*,
double*),
const int& n,
double* A,
const int& lda,
double* B,
const int& ldb,
int* sdim,
double* ALPHAR,
double* ALPHAI,
double* BETA,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
double* WORK,
const int& lwork,
int* bwork,
int* info )
const;
1667 void GESVD(
const char& JOBU,
const char& JOBVT,
const int& m,
const int& n,
double* A,
const int& lda,
double* S,
double* U,
const int& ldu,
double* V,
const int& ldv,
double* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1670 void ORMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
double* A,
const int& lda,
const double* TAU,
double* C,
const int& ldc,
double* WORK,
const int& lwork,
int* info)
const;
1671 void ORM2R(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
const double A[],
const int& lda,
const double TAU[],
double C[],
const int& ldc,
double WORK[],
int*
const info)
const;
1672 void UNMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
double* A,
const int& lda,
const double* TAU,
double* C,
const int& ldc,
double* WORK,
const int& lwork,
int* info)
const;
1673 void UNM2R(
const char& SIDE,
const char&
TRANS,
const int& M,
const int& N,
const int& K,
const double A[],
const int& LDA,
const double TAU[],
double C[],
const int& LDC,
double WORK[],
int*
const INFO)
const;
1674 void ORGQR(
const int& m,
const int& n,
const int& k,
double* A,
const int& lda,
const double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1675 void UNGQR(
const int& m,
const int& n,
const int& k,
double* A,
const int& lda,
const double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1676 void ORGHR(
const int& n,
const int& ilo,
const int& ihi,
double* A,
const int& lda,
const double* TAU,
double* WORK,
const int& lwork,
int* info)
const;
1677 void ORMHR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& ilo,
const int& ihi,
const double* A,
const int& lda,
const double* TAU,
double* C,
const int& ldc,
double* WORK,
const int& lwork,
int* info)
const;
1680 void TREVC(
const char& SIDE,
const char& HOWMNY,
int* select,
const int& n,
const double* T,
const int& ldt,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
const int& mm,
int* m,
double* WORK,
int* info)
const;
1681 void TREVC(
const char& SIDE,
const int& n,
const double* T,
const int& ldt,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
const int& mm,
int* m,
double* WORK,
double* RWORK,
int* info)
const;
1684 TEUCHOS_DEPRECATED
void TREXC(
const char& COMPQ,
const int& n,
double* T,
const int& ldt,
double* Q,
const int& ldq,
const int& ifst,
const int& ilst,
double* WORK,
int* info)
const;
1685 void TREXC(
const char& COMPQ,
const int& n,
double* T,
const int& ldt,
double* Q,
const int& ldq,
int* ifst,
int* ilst,
double* WORK,
int* info)
const;
1687 void TGEVC(
const char& SIDE,
const char& HOWMNY,
const int* SELECT,
const int& n,
double* S,
const int& lds,
double* P,
const int& ldp,
double* VL,
const int& ldvl,
double* VR,
const int& ldvr,
const int& mm,
int* M,
double* WORK,
int* info)
const;
1690 void LARTG(
const double& f,
const double& g,
double* c,
double* s,
double* r )
const;
1691 void LARFG(
const int& n,
double* alpha,
double* x,
const int& incx,
double* tau )
const;
1696 TEUCHOS_DEPRECATED
void GEBAL(
const char& JOBZ,
const int& n,
double* A,
const int& lda,
const int& ilo,
const int& ihi,
double* scale,
int* info)
const;
1697 void GEBAL(
const char& JOBZ,
const int& n,
double* A,
const int& lda,
int* ilo,
int* ihi,
double* scale,
int* info)
const;
1699 void GEBAK(
const char& JOBZ,
const char& SIDE,
const int& n,
const int& ilo,
const int& ihi,
const double* scale,
const int& m,
double* V,
const int& ldv,
int* info)
const;
1702 double LARND(
const int& idist,
int* seed )
const;
1703 void LARNV(
const int& idist,
int* seed,
const int& n,
double* v )
const;
1706 double LAMCH(
const char& CMACH)
const;
1707 int ILAENV(
const int& ispec,
const std::string& NAME,
const std::string& OPTS,
const int& N1 = -1,
const int& N2 = -1,
const int& N3 = -1,
const int& N4 = -1 )
const;
1710 double LAPY2(
const double& x,
const double& y)
const;
1716 #ifdef HAVE_TEUCHOS_COMPLEX
1721 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, std::complex<float> >
1724 inline LAPACK(
void) {}
1725 inline LAPACK(
const LAPACK<
int, std::complex<float> >& lapack) {}
1726 inline virtual ~LAPACK(
void) {}
1729 void PTTRF(
const int& n, std::complex<float>* d, std::complex<float>* e,
int* info)
const;
1730 void PTTRS(
const int& n,
const int& nrhs,
const std::complex<float>* d,
const std::complex<float>* e, std::complex<float>* B,
const int& ldb,
int* info)
const;
1731 void POTRF(
const char& UPLO,
const int& n, std::complex<float>* A,
const int& lda,
int* info)
const;
1732 void POTRS(
const char& UPLO,
const int& n,
const int& nrhs,
const std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
int* info)
const;
1733 void POTRI(
const char& UPLO,
const int& n, std::complex<float>* A,
const int& lda,
int* info)
const;
1734 void POCON(
const char& UPLO,
const int& n,
const std::complex<float>* A,
const int& lda,
const float& anorm,
float* rcond, std::complex<float>* WORK,
float* rwork,
int* info)
const;
1735 void POSV(
const char& UPLO,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
int* info)
const;
1736 void POEQU(
const int& n,
const std::complex<float>* A,
const int& lda,
float* S,
float* scond,
float* amax,
int* info)
const;
1737 void PORFS(
const char& UPLO,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda,
const std::complex<float>* AF,
const int& ldaf,
const std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1740 TEUCHOS_DEPRECATED
void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* AF,
const int& ldaf,
const char& EQUED,
float* S, std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1741 void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* AF,
const int& ldaf,
char* EQUED,
float* S, std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1744 void GELS(
const char&
TRANS,
const int& m,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1745 void GELSS(
const int& m,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
float* S,
const float& rcond,
int* rank, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1746 void GEQRF(
const int& m,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* TAU, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1747 void GEQR2(
const int& m,
const int& n, std::complex<float> A[],
const int& lda, std::complex<float> TAU[], std::complex<float> WORK[],
int*
const info)
const;
1748 void UNGQR(
const int& m,
const int& n,
const int& k, std::complex<float>* A,
const int& lda,
const std::complex<float>* TAU, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1749 void UNMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k, std::complex<float>* A,
const int& lda,
const std::complex<float>* TAU, std::complex<float>* C,
const int& ldc, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1750 void UNM2R(
const char& SIDE,
const char&
TRANS,
const int& M,
const int& N,
const int& K,
const std::complex<float> A[],
const int& LDA,
const std::complex<float> TAU[], std::complex<float> C[],
const int& LDC, std::complex<float> WORK[],
int*
const INFO)
const;
1751 void GETRF(
const int& m,
const int& n, std::complex<float>* A,
const int& lda,
int* IPIV,
int* info)
const;
1752 void GETRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<float>* A,
const int& lda,
const int* IPIV, std::complex<float>* B,
const int& ldb,
int* info)
const;
1753 void LASCL(
const char& TYPE,
const int& kl,
const int& ku,
const float& cfrom,
const float& cto,
const int& m,
const int& n, std::complex<float>* A,
const int& lda,
int* info)
const;
1756 GEQP3 (
const int& m,
1758 std::complex<float>* A,
1761 std::complex<float>* TAU,
1762 std::complex<float>* WORK,
1767 void LASWP (
const int& N,
1768 std::complex<float> A[],
1773 const int& INCX)
const;
1775 void GBTRF(
const int& m,
const int& n,
const int& kl,
const int& ku, std::complex<float>* A,
const int& lda,
int* IPIV,
int* info)
const;
1776 void GBTRS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const std::complex<float>* A,
const int& lda,
const int* IPIV, std::complex<float>* B,
const int& ldb,
int* info)
const;
1777 void GTTRF(
const int& n, std::complex<float>* dl, std::complex<float>* d, std::complex<float>* du, std::complex<float>* du2,
int* IPIV,
int* info)
const;
1778 void GTTRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<float>* dl,
const std::complex<float>* d,
const std::complex<float>* du,
const std::complex<float>* du2,
const int* IPIV, std::complex<float>* B,
const int& ldb,
int* info)
const;
1779 void GETRI(
const int& n, std::complex<float>* A,
const int& lda,
const int* IPIV, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1780 void LATRS (
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const char& NORMIN,
const int& N, std::complex<float>* A,
const int& LDA, std::complex<float>* X,
float* SCALE,
float* CNORM,
int* INFO)
const;
1781 void GECON(
const char& NORM,
const int& n,
const std::complex<float>* A,
const int& lda,
const float& anorm,
float* rcond, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1782 void GBCON(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const std::complex<float>* A,
const int& lda,
int* IPIV,
const float& anorm,
float* rcond, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1783 float LANGB(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const std::complex<float>* A,
const int& lda,
float* WORK)
const;
1784 void GESV(
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda,
int* IPIV, std::complex<float>* B,
const int& ldb,
int* info)
const;
1785 void GEEQU(
const int& m,
const int& n,
const std::complex<float>* A,
const int& lda,
float* R,
float* C,
float* rowcond,
float* colcond,
float* amax,
int* info)
const;
1786 void GERFS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<float>* A,
const int& lda,
const std::complex<float>* AF,
const int& ldaf,
const int* IPIV,
const std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1787 void GBEQU(
const int& m,
const int& n,
const int& kl,
const int& ku,
const std::complex<float>* A,
const int& lda,
float* R,
float* C,
float* rowcond,
float* colcond,
float* amax,
int* info)
const;
1788 void GBRFS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const std::complex<float>* A,
const int& lda,
const std::complex<float>* AF,
const int& ldaf,
const int* IPIV,
const std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1791 TEUCHOS_DEPRECATED
void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* AF,
const int& ldaf,
int* IPIV,
const char& EQUED,
float* R,
float* C, std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1792 void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs, std::complex<float>* A,
const int& lda, std::complex<float>* AF,
const int& ldaf,
int* IPIV,
char* EQUED,
float* R,
float* C, std::complex<float>* B,
const int& ldb, std::complex<float>* X,
const int& ldx,
float* rcond,
float* FERR,
float* BERR, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1794 void GEHRD(
const int& n,
const int& ilo,
const int& ihi, std::complex<float>* A,
const int& lda, std::complex<float>* TAU, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1795 void TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const int& n,
const int& nrhs,
const std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
int* info)
const;
1796 void TRTRI(
const char& UPLO,
const char& DIAG,
const int& n,
const std::complex<float>* A,
const int& lda,
int* info)
const;
1799 void STEQR(
const char& COMPZ,
const int& n,
float* D,
float* E, std::complex<float>* Z,
const int& ldz,
float* WORK,
int* info)
const;
1800 void HEEV(
const char& JOBZ,
const char& UPLO,
const int& n, std::complex<float>* A,
const int& lda,
float* W, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1801 void HEGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
float* W, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1804 void HSEQR(
const char& JOB,
const char& COMPZ,
const int& n,
const int& ilo,
const int& ihi, std::complex<float>* H,
const int& ldh, std::complex<float>* W, std::complex<float>* Z,
const int& ldz, std::complex<float>* WORK,
const int& lwork,
int* info)
const;
1805 void GEES(
const char& JOBVS,
const char& SORT,
int (*ptr2func)(std::complex<float>*),
const int& n, std::complex<float>* A,
const int& lda,
int* sdim, std::complex<float>* W, std::complex<float>* VS,
const int& ldvs, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* BWORK,
int* info)
const;
1806 void GEES(
const char& JOBVS,
const int& n, std::complex<float>* A,
const int& lda,
int* sdim,
float* WR,
float* WI, std::complex<float>* VS,
const int& ldvs, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* BWORK,
int* info)
const;
1808 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* W, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1809 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<float>* A,
const int& lda,
float* WR,
float* WI, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1811 void GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* W, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* SCALE,
float* abnrm,
float* RCONDE,
float* RCONDV, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1813 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb, std::complex<float>* ALPHA, std::complex<float>* BETA, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* lscale,
float* rscale,
float* abnrm,
float* bbnrm,
float* RCONDE,
float* RCONDV, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* IWORK,
int* BWORK,
int* info)
const;
1814 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<float>* A,
const int& lda, std::complex<float>* B,
const int& ldb,
float* ALPHAR,
float* ALPHAI, std::complex<float>* BETA, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr,
int* ilo,
int* ihi,
float* lscale,
float* rscale,
float* abnrm,
float* bbnrm,
float* RCONDE,
float* RCONDV, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* IWORK,
int* BWORK,
int* info)
const;
1815 void GGEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<float> *A,
const int& lda, std::complex<float> *B,
const int& ldb, std::complex<float>* ALPHA, std::complex<float>* BETA, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr, std::complex<float> *WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1818 void GESVD(
const char& JOBU,
const char& JOBVT,
const int& m,
const int& n, std::complex<float>* A,
const int& lda,
float* S, std::complex<float>* U,
const int& ldu, std::complex<float>* V,
const int& ldv, std::complex<float>* WORK,
const int& lwork,
float* RWORK,
int* info)
const;
1821 void TREVC(
const char& SIDE,
const char& HOWMNY,
int* select,
const int& n,
const std::complex<float>* T,
const int& ldt, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr,
const int& mm,
int* m, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1822 void TREVC(
const char& SIDE,
const int& n,
const std::complex<float>* T,
const int& ldt, std::complex<float>* VL,
const int& ldvl, std::complex<float>* VR,
const int& ldvr,
const int& mm,
int* m, std::complex<float>* WORK,
float* RWORK,
int* info)
const;
1825 TEUCHOS_DEPRECATED
void TREXC(
const char& COMPQ,
const int& n, std::complex<float>* T,
const int& ldt, std::complex<float>* Q,
const int& ldq,
const int& ifst,
const int& ilst, std::complex<float>* WORK,
int* info)
const;
1826 void TREXC(
const char& COMPQ,
const int& n, std::complex<float>* T,
const int& ldt, std::complex<float>* Q,
const int& ldq,
int* ifst,
int* ilst, std::complex<float>* WORK,
int* info)
const;
1829 void LARTG(
const std::complex<float> f,
const std::complex<float> g,
float* c, std::complex<float>* s, std::complex<float>* r )
const;
1830 void LARFG(
const int& n, std::complex<float>* alpha, std::complex<float>* x,
const int& incx, std::complex<float>* tau )
const;
1835 TEUCHOS_DEPRECATED
void GEBAL(
const char& JOBZ,
const int& n, std::complex<float>* A,
const int& lda,
const int& ilo,
const int& ihi,
float* scale,
int* info)
const;
1836 void GEBAL(
const char& JOBZ,
const int& n, std::complex<float>* A,
const int& lda,
int* ilo,
int* ihi,
float* scale,
int* info)
const;
1838 void GEBAK(
const char& JOBZ,
const char& SIDE,
const int& n,
const int& ilo,
const int& ihi,
const float* scale,
const int& m, std::complex<float>* V,
const int& ldv,
int* info)
const;
1841 std::complex<float> LARND(
const int& idist,
int* seed )
const;
1842 void LARNV(
const int& idist,
int* seed,
const int& n, std::complex<float>* v )
const;
1845 int ILAENV(
const int& ispec,
const std::string& NAME,
const std::string& OPTS,
const int& N1 = -1,
const int& N2 = -1,
const int& N3 = -1,
const int& N4 = -1 )
const;
1854 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, std::complex<double> >
1857 inline LAPACK(
void) {}
1858 inline LAPACK(
const LAPACK<
int, std::complex<double> >& lapack) {}
1859 inline virtual ~LAPACK(
void) {}
1862 void PTTRF(
const int& n, std::complex<double>* d, std::complex<double>* e,
int* info)
const;
1863 void PTTRS(
const int& n,
const int& nrhs,
const std::complex<double>* d,
const std::complex<double>* e, std::complex<double>* B,
const int& ldb,
int* info)
const;
1864 void POTRF(
const char& UPLO,
const int& n, std::complex<double>* A,
const int& lda,
int* info)
const;
1865 void POTRS(
const char& UPLO,
const int& n,
const int& nrhs,
const std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
int* info)
const;
1866 void POTRI(
const char& UPLO,
const int& n, std::complex<double>* A,
const int& lda,
int* info)
const;
1867 void POCON(
const char& UPLO,
const int& n,
const std::complex<double>* A,
const int& lda,
const double& anorm,
double* rcond, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1868 void POSV(
const char& UPLO,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
int* info)
const;
1869 void POEQU(
const int& n,
const std::complex<double>* A,
const int& lda,
double* S,
double* scond,
double* amax,
int* info)
const;
1870 void PORFS(
const char& UPLO,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda,
const std::complex<double>* AF,
const int& ldaf,
const std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1873 TEUCHOS_DEPRECATED
void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* AF,
const int& ldaf,
const char& EQUED,
double* S, std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1874 void POSVX(
const char& FACT,
const char& UPLO,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* AF,
const int& ldaf,
char* EQUED,
double* S, std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1877 void GELS(
const char&
TRANS,
const int& m,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1878 void GELSS(
const int& m,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
double* S,
const double& rcond,
int* rank, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1879 void GEQRF(
const int& m,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* TAU, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1880 void GEQR2(
const int& m,
const int& n, std::complex<double> A[],
const int& lda, std::complex<double> TAU[], std::complex<double> WORK[],
int*
const info)
const;
1881 void UNGQR(
const int& m,
const int& n,
const int& k, std::complex<double>* A,
const int& lda,
const std::complex<double>* TAU, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1882 void UNMQR(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k, std::complex<double>* A,
const int& lda,
const std::complex<double>* TAU, std::complex<double>* C,
const int& ldc, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1883 void UNM2R(
const char& SIDE,
const char&
TRANS,
const int& M,
const int& N,
const int& K,
const std::complex<double> A[],
const int& LDA,
const std::complex<double> TAU[], std::complex<double> C[],
const int& LDC, std::complex<double> WORK[],
int*
const INFO)
const;
1885 void GETRF(
const int& m,
const int& n, std::complex<double>* A,
const int& lda,
int* IPIV,
int* info)
const;
1886 void GETRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<double>* A,
const int& lda,
const int* IPIV, std::complex<double>* B,
const int& ldb,
int* info)
const;
1887 void LASCL(
const char& TYPE,
const int& kl,
const int& ku,
const double& cfrom,
const double& cto,
const int& m,
const int& n, std::complex<double>* A,
const int& lda,
int* info)
const;
1890 GEQP3 (
const int& m,
1892 std::complex<double>* A,
1895 std::complex<double>* TAU,
1896 std::complex<double>* WORK,
1901 void LASWP (
const int& N,
1902 std::complex<double> A[],
1907 const int& INCX)
const;
1909 void GBTRF(
const int& m,
const int& n,
const int& kl,
const int& ku, std::complex<double>* A,
const int& lda,
int* IPIV,
int* info)
const;
1910 void GBTRS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const std::complex<double>* A,
const int& lda,
const int* IPIV, std::complex<double>* B,
const int& ldb,
int* info)
const;
1911 void GTTRF(
const int& n, std::complex<double>* dl, std::complex<double>* d, std::complex<double>* du, std::complex<double>* du2,
int* IPIV,
int* info)
const;
1912 void GTTRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<double>* dl,
const std::complex<double>* d,
const std::complex<double>* du,
const std::complex<double>* du2,
const int* IPIV, std::complex<double>* B,
const int& ldb,
int* info)
const;
1913 void GETRI(
const int& n, std::complex<double>* A,
const int& lda,
const int* IPIV, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1914 void LATRS (
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const char& NORMIN,
const int& N, std::complex<double>* A,
const int& LDA, std::complex<double>* X,
double* SCALE,
double* CNORM,
int* INFO)
const;
1915 void GECON(
const char& NORM,
const int& n,
const std::complex<double>* A,
const int& lda,
const double& anorm,
double* rcond, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1916 void GBCON(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const std::complex<double>* A,
const int& lda,
int* IPIV,
const double& anorm,
double* rcond, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1917 double LANGB(
const char& NORM,
const int& n,
const int& kl,
const int& ku,
const std::complex<double>* A,
const int& lda,
double* WORK)
const;
1918 void GESV(
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda,
int* IPIV, std::complex<double>* B,
const int& ldb,
int* info)
const;
1919 void GEEQU(
const int& m,
const int& n,
const std::complex<double>* A,
const int& lda,
double* R,
double* C,
double* rowcond,
double* colcond,
double* amax,
int* info)
const;
1920 void GERFS(
const char&
TRANS,
const int& n,
const int& nrhs,
const std::complex<double>* A,
const int& lda,
const std::complex<double>* AF,
const int& ldaf,
const int* IPIV,
const std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1921 void GBEQU(
const int& m,
const int& n,
const int& kl,
const int& ku,
const std::complex<double>* A,
const int& lda,
double* R,
double* C,
double* rowcond,
double* colcond,
double* amax,
int* info)
const;
1922 void GBRFS(
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const std::complex<double>* A,
const int& lda,
const std::complex<double>* AF,
const int& ldaf,
const int* IPIV,
const std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1925 TEUCHOS_DEPRECATED
void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* AF,
const int& ldaf,
int* IPIV,
const char& EQUED,
double* R,
double* C, std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1926 void GESVX(
const char& FACT,
const char&
TRANS,
const int& n,
const int& nrhs, std::complex<double>* A,
const int& lda, std::complex<double>* AF,
const int& ldaf,
int* IPIV,
char* EQUED,
double* R,
double* C, std::complex<double>* B,
const int& ldb, std::complex<double>* X,
const int& ldx,
double* rcond,
double* FERR,
double* BERR, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1928 void GEHRD(
const int& n,
const int& ilo,
const int& ihi, std::complex<double>* A,
const int& lda, std::complex<double>* TAU, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1929 void TRTRS(
const char& UPLO,
const char&
TRANS,
const char& DIAG,
const int& n,
const int& nrhs,
const std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
int* info)
const;
1930 void TRTRI(
const char& UPLO,
const char& DIAG,
const int& n,
const std::complex<double>* A,
const int& lda,
int* info)
const;
1933 void STEQR(
const char& COMPZ,
const int& n,
double* D,
double* E, std::complex<double>* Z,
const int& ldz,
double* WORK,
int* info)
const;
1934 void HEEV(
const char& JOBZ,
const char& UPLO,
const int& n, std::complex<double>* A,
const int& lda,
double* W, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1935 void HEGV(
const int& itype,
const char& JOBZ,
const char& UPLO,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
double* W, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1938 void HSEQR(
const char& JOB,
const char& COMPZ,
const int& n,
const int& ilo,
const int& ihi, std::complex<double>* H,
const int& ldh, std::complex<double>* W, std::complex<double>* Z,
const int& ldz, std::complex<double>* WORK,
const int& lwork,
int* info)
const;
1939 void GEES(
const char& JOBVS,
const char& SORT,
int (*ptr2func)(std::complex<double>*),
const int& n, std::complex<double>* A,
const int& lda,
int* sdim, std::complex<double>* W, std::complex<double>* VS,
const int& ldvs, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* BWORK,
int* info)
const;
1940 void GEES(
const char& JOBVS,
const int& n, std::complex<double>* A,
const int& lda,
int* sdim,
double* WR,
double* WI, std::complex<double>* VS,
const int& ldvs, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* BWORK,
int* info)
const;
1942 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* W, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1943 void GEEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<double>* A,
const int& lda,
double* WR,
double* WI, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1945 void GEEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* W, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* SCALE,
double* abnrm,
double* RCONDE,
double* RCONDV, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1946 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb, std::complex<double>* ALPHA, std::complex<double>* BETA, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* lscale,
double* rscale,
double* abnrm,
double* bbnrm,
double* RCONDE,
double* RCONDV, std::complex<double>* work,
const int& lwork,
double* RWORK,
int* IWORK,
int* BWORK,
int* info)
const;
1947 void GGEVX(
const char& BALANC,
const char& JOBVL,
const char& JOBVR,
const char& SENSE,
const int& n, std::complex<double>* A,
const int& lda, std::complex<double>* B,
const int& ldb,
double* ALPHAR,
double* ALPHAI, std::complex<double>* BETA, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr,
int* ilo,
int* ihi,
double* lscale,
double* rscale,
double* abnrm,
double* bbnrm,
double* RCONDE,
double* RCONDV, std::complex<double>* work,
const int& lwork,
double* RWORK,
int* IWORK,
int* BWORK,
int* info)
const;
1948 void GGEV(
const char& JOBVL,
const char& JOBVR,
const int& n, std::complex<double> *A,
const int& lda, std::complex<double> *B,
const int& ldb, std::complex<double>* ALPHA, std::complex<double>* BETA, std::complex<double>* VL,
const int& ldvl, std::complex<double>*VR,
const int& ldvr, std::complex<double> *WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1951 void GESVD(
const char& JOBU,
const char& JOBVT,
const int& m,
const int& n, std::complex<double>* A,
const int& lda,
double* S, std::complex<double>* U,
const int& ldu, std::complex<double>* V,
const int& ldv, std::complex<double>* WORK,
const int& lwork,
double* RWORK,
int* info)
const;
1954 void TREVC(
const char& SIDE,
const char& HOWMNY,
int* select,
const int& n,
const std::complex<double>* T,
const int& ldt, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr,
const int& mm,
int* m, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1955 void TREVC(
const char& SIDE,
const int& n,
const std::complex<double>* T,
const int& ldt, std::complex<double>* VL,
const int& ldvl, std::complex<double>* VR,
const int& ldvr,
const int& mm,
int* m, std::complex<double>* WORK,
double* RWORK,
int* info)
const;
1958 TEUCHOS_DEPRECATED
void TREXC(
const char& COMPQ,
const int& n, std::complex<double>* T,
const int& ldt, std::complex<double>* Q,
const int& ldq,
const int& ifst,
const int& ilst, std::complex<double>* WORK,
int* info)
const;
1959 void TREXC(
const char& COMPQ,
const int& n, std::complex<double>* T,
const int& ldt, std::complex<double>* Q,
const int& ldq,
int* ifst,
int* ilst, std::complex<double>* WORK,
int* info)
const;
1962 void LARTG(
const std::complex<double> f,
const std::complex<double> g,
double* c, std::complex<double>* s, std::complex<double>* r )
const;
1963 void LARFG(
const int& n, std::complex<double>* alpha, std::complex<double>* x,
const int& incx, std::complex<double>* tau )
const;
1968 TEUCHOS_DEPRECATED
void GEBAL(
const char& JOBZ,
const int& n, std::complex<double>* A,
const int& lda,
const int& ilo,
const int& ihi,
double* scale,
int* info)
const;
1969 void GEBAL(
const char& JOBZ,
const int& n, std::complex<double>* A,
const int& lda,
int* ilo,
int* ihi,
double* scale,
int* info)
const;
1971 void GEBAK(
const char& JOBZ,
const char& SIDE,
const int& n,
const int& ilo,
const int& ihi,
const double* scale,
const int& m, std::complex<double>* V,
const int& ldv,
int* info)
const;
1974 std::complex<double> LARND(
const int& idist,
int* seed )
const;
1975 void LARNV(
const int& idist,
int* seed,
const int& n, std::complex<double>* v )
const;
1978 int ILAENV(
const int& ispec,
const std::string& NAME,
const std::string& OPTS,
const int& N1 = -1,
const int& N2 = -1,
const int& N3 = -1,
const int& N4 = -1 )
const;
1984 #endif // HAVE_TEUCHOS_COMPLEX
1986 #ifdef HAVE_TEUCHOSCORE_QUADMATH
1993 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, __float128>
1996 inline LAPACK(
void) {}
1997 inline LAPACK(
const LAPACK<int, __float128>& lapack) {}
1998 inline virtual ~LAPACK(
void) {}
2000 void GEQRF(
const int& m,
const int& n, __float128* A,
const int& lda, __float128* TAU, __float128* WORK,
const int& lwork,
int* info)
const;
2001 void GEQR2(
const int& m,
const int& n, __float128 A[],
const int& lda, __float128 TAU[], __float128 WORK[],
int*
const info)
const;
2002 void GETRF(
const int& m,
const int& n, __float128* A,
const int& lda,
int* IPIV,
int* info)
const;
2003 void GETRS(
const char&
TRANS,
const int& n,
const int& nrhs,
const __float128* A,
const int& lda,
const int* IPIV, __float128* B,
const int& ldb,
int* info)
const;
2004 void GETRI(
const int& n, __float128* A,
const int& lda,
const int* IPIV, __float128* WORK,
const int& lwork,
int* info)
const;
2005 void LASWP (
const int& N, __float128 A[],
const int& LDA,
const int& K1,
const int& K2,
const int IPIV[],
const int& INCX)
const;
2007 void ORM2R(
const char& SIDE,
const char&
TRANS,
const int& m,
const int& n,
const int& k,
const __float128 A[],
const int& lda,
const __float128 TAU[], __float128 C[],
const int& ldc, __float128 WORK[],
int*
const info)
const;
2008 void ORGQR(
const int& m,
const int& n,
const int& k, __float128* A,
const int& lda,
const __float128* TAU, __float128* WORK,
const int& lwork,
int* info)
const;
2009 void UNGQR(
const int& m,
const int& n,
const int& k, __float128* A,
const int& lda,
const __float128* TAU, __float128* WORK,
const int& lwork,
int* info)
const;
2011 void LARFG(
const int& n, __float128* alpha, __float128* x,
const int& incx, __float128* tau )
const;
2013 __float128 LAPY2 (
const __float128 x,
const __float128 y)
const;
2014 void LASCL (
const char& TYPE,
const int& kl,
const int& ku,
const __float128 cfrom,
const __float128 cto,
const int& m,
const int& n, __float128* A,
const int& lda,
int* info)
const;
2016 void GBTRF (
const int& m,
const int& n,
const int& kl,
const int& ku, __float128* A,
const int& lda,
int* IPIV,
int* info)
const;
2017 void GBTRS (
const char&
TRANS,
const int& n,
const int& kl,
const int& ku,
const int& nrhs,
const __float128* A,
const int& lda,
const int* IPIV, __float128* B,
const int& ldb,
int* info)
const;
2022 #endif // HAVE_TEUCHOSCORE_QUADMATH
2024 #endif // DOXYGEN_SHOULD_SKIP_THIS
2028 #endif // _TEUCHOS_LAPACK_HPP_