42 #ifndef _TEUCHOS_LAPACK_UQ_PCE_HPP_
43 #define _TEUCHOS_LAPACK_UQ_PCE_HPP_
58 template<
typename OrdinalType,
typename Storage>
59 class LAPACK<OrdinalType, Sacado::UQ::PCE<Storage> >
84 { throw_error(
"PTTRF"); }
88 { throw_error(
"PTTRS"); }
91 void POTRF(
const char UPLO,
const OrdinalType n,
ScalarType*
A,
const OrdinalType lda, OrdinalType* info)
const
92 { throw_error(
"POTRF"); }
95 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
96 { throw_error(
"POTRS"); }
99 void POTRI(
const char UPLO,
const OrdinalType n,
ScalarType*
A,
const OrdinalType lda, OrdinalType* info)
const
100 { throw_error(
"POTRI"); }
105 { throw_error(
"POCON"); }
108 void POSV(
const char UPLO,
const OrdinalType n,
const OrdinalType nrhs,
ScalarType*
A,
const OrdinalType lda,
ScalarType*
B,
const OrdinalType ldb, OrdinalType* info)
const
109 { throw_error(
"POSV"); }
113 { throw_error(
"POEQU"); }
116 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
117 { throw_error(
"PORFS"); }
120 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
121 { throw_error(
"POSVX"); }
128 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
129 { throw_error(
"GELS"); }
164 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
165 { throw_error(
"GELSS"); }
168 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
169 { throw_error(
"GELSS"); }
172 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
173 { throw_error(
"GGLSE"); }
177 { throw_error(
"GEQRF"); }
180 void GETRF(
const OrdinalType m,
const OrdinalType n,
ScalarType*
A,
const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info)
const
181 { throw_error(
"GETRF"); }
184 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
185 { throw_error(
"GETRS"); }
188 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
189 { throw_error(
"LASCL"); }
195 const OrdinalType lda,
199 const OrdinalType lwork,
201 OrdinalType* info )
const
202 { throw_error(
"GEQP3"); }
208 const OrdinalType LDA,
209 const OrdinalType K1,
210 const OrdinalType K2,
211 const OrdinalType IPIV[],
212 const OrdinalType INCX)
const
213 { throw_error(
"LASWP"); }
216 void GBTRF(
const OrdinalType m,
const OrdinalType n,
const OrdinalType kl,
const OrdinalType ku,
ScalarType*
A,
const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info)
const
217 { throw_error(
"GBTRF"); }
220 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
221 { throw_error(
"GBTRS"); }
225 { throw_error(
"GTTRF"); }
228 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
229 { throw_error(
"GTTRS"); }
232 void GETRI(
const OrdinalType n,
ScalarType*
A,
const OrdinalType lda,
const OrdinalType* IPIV,
ScalarType* WORK,
const OrdinalType lwork, OrdinalType* info)
const
233 { throw_error(
"GETRI"); }
246 const OrdinalType LDA,
250 OrdinalType* INFO)
const
251 { throw_error(
"LATRS"); }
255 { throw_error(
"GECON"); }
258 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
259 { throw_error(
"GBCON"); }
266 void GESV(
const OrdinalType n,
const OrdinalType nrhs,
ScalarType*
A,
const OrdinalType lda, OrdinalType* IPIV,
ScalarType*
B,
const OrdinalType ldb, OrdinalType* info)
const
267 { throw_error(
"GESV"); }
271 { throw_error(
"GEEQU"); }
274 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,
MagnitudeType* FERR,
MagnitudeType* BERR,
ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info)
const
275 { throw_error(
"GERFS"); }
278 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
279 { throw_error(
"GBEQU"); }
282 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
283 { throw_error(
"GBRFS"); }
286 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
287 { throw_error(
"GESVX"); }
293 { throw_error(
"SYTRD"); }
296 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
297 { throw_error(
"GEHRD"); }
300 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
301 { throw_error(
"TRTRS"); }
304 void TRTRI(
const char UPLO,
const char DIAG,
const OrdinalType n,
const ScalarType*
A,
const OrdinalType lda, OrdinalType* info)
const
305 { throw_error(
"TRTRI"); }
314 { throw_error(
"SPEV"); }
320 { throw_error(
"SYEV"); }
325 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
326 { throw_error(
"SYGV"); }
332 { throw_error(
"HEEV"); }
337 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
338 { throw_error(
"HEGV"); }
342 { throw_error(
"STEQR"); }
347 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
349 { throw_error(
"HSEQR"); }
354 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
355 { throw_error(
"GEES"); }
360 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
361 { throw_error(
"GEES"); }
366 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
367 { throw_error(
"GEES"); }
374 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 { throw_error(
"GEEV"); }
381 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
382 { throw_error(
"GEEVX"); }
388 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
389 { throw_error(
"GGEVX"); }
394 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
395 { throw_error(
"GGEV"); }
401 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
402 { throw_error(
"TRSEN"); }
408 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
409 { throw_error(
"TGSEN"); }
415 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
416 { throw_error(
"GGES"); }
423 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
425 { throw_error(
"GESVD"); }
441 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
442 { throw_error(
"ORMQR"); }
452 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
453 { throw_error(
"UNMQR"); }
464 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
465 { throw_error(
"ORGQR"); }
475 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
476 { throw_error(
"UNGQR"); }
481 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 { throw_error(
"ORGHR"); }
487 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
488 { throw_error(
"ORMHR"); }
496 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
497 { throw_error(
"TREVC"); }
502 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
503 { throw_error(
"TREVC"); }
508 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
509 { throw_error(
"TREXC"); }
514 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
515 { throw_error(
"TGEVC"); }
525 { throw_error(
"LARTG"); }
529 { throw_error(
"LARFG"); }
537 void GEBAL(
const char JOBZ,
const OrdinalType n,
ScalarType*
A,
const OrdinalType lda, OrdinalType ilo, OrdinalType ihi,
MagnitudeType* scale, OrdinalType* info)
const
538 { throw_error(
"GEBAL"); }
541 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
542 { throw_error(
"GEBAK"); }
548 ScalarType
LARND(
const OrdinalType idist, OrdinalType* seed )
const
550 { throw_error(
"LARND");
return ScalarType(0); }
553 void LARNV(
const OrdinalType idist, OrdinalType* seed,
const OrdinalType n,
ScalarType* v )
const
554 { throw_error(
"LARNV"); }
563 { throw_error(
"LAMCH");
return ScalarType(0); }
569 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
570 { throw_error(
"ILAENV");
return OrdinalType(0); }
579 { throw_error(
"LAPY2");
return ScalarType(0); }
586 true, std::logic_error,
587 func <<
": Not implemented for Sacado::UQ::PCE scalar type!");
594 #endif // _TEUCHOS_LAPACK__UQ_PCE_HPP_
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
void LASWP(const OrdinalType N, ScalarType A[], const OrdinalType LDA, const OrdinalType K1, const OrdinalType K2, const OrdinalType IPIV[], const OrdinalType INCX) const
Apply a series of row interchanges to the matrix A.
void PTTRF(const OrdinalType n, ScalarType *d, ScalarType *e, OrdinalType *info) const
Computes the L*D*L' factorization of a Hermitian/symmetric positive definite tridiagonal matrix A...
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
Uses the Cholesky factorization to compute the solution to a real system of linear equations A*X=B...
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
Computes all the eigenvalues and, optionally, eigenvectors of a Hermitian n by n matrix A...
void STEQR(const char COMPZ, const OrdinalType n, ScalarType *D, ScalarType *E, ScalarType *Z, const OrdinalType ldz, ScalarType *WORK, OrdinalType *info) const
Computes the eigenvalues and, optionally, eigenvectors of a symmetric tridiagonal n by n matrix A usi...
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
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
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix pencil {A...
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
void GBTRF(const OrdinalType m, const OrdinalType n, const OrdinalType kl, const OrdinalType ku, ScalarType *A, const OrdinalType lda, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a general banded m by n matrix A using partial pivoting with row inte...
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
Apply Householder reflectors (complex case).
void GETRI(const OrdinalType n, ScalarType *A, const OrdinalType lda, const OrdinalType *IPIV, ScalarType *WORK, const OrdinalType lwork, OrdinalType *info) const
Computes the inverse of a matrix A using the LU factorization computed by GETRF.
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
Legacy GELSS interface for real-valued ScalarType.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual ~LAPACK(void)
Destructor.
void throw_error(const char *func) const
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
void POTRI(const char UPLO, const OrdinalType n, ScalarType *A, const OrdinalType lda, OrdinalType *info) const
Computes the inverse of a real symmetric positive definite matrix A using the Cholesky factorization ...
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
void GEQP3(const OrdinalType m, const OrdinalType n, ScalarType *A, const OrdinalType lda, OrdinalType *jpvt, ScalarType *TAU, ScalarType *WORK, const OrdinalType lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes a QR factorization with column pivoting of a matrix A: A*P = Q*R using Level 3 BLAS...
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
Use the SVD to solve a possibly rank-deficient linear least-squares problem.
void LARNV(const OrdinalType idist, OrdinalType *seed, const OrdinalType n, ScalarType *v) const
Returns a vector of random numbers from a chosen distribution.
void LARFG(const OrdinalType n, ScalarType *alpha, ScalarType *x, const OrdinalType incx, ScalarType *tau) const
Generates an elementary reflector of order n that zeros out the last n-1 components of the input vect...
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
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
Forms the left or right eigenvectors of a general matrix that has been balanced by GEBAL by backward ...
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
Estimates the reciprocal of the condition number (1-norm) of a real symmetric positive definite matri...
Sacado::UQ::PCE< Storage > ScalarType
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
Overwrites the general real m by n matrix C with the product of C and Q, which is a product of ihi-il...
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
Computes all the eigenvalues and, optionally, eigenvectors of a generalized Hermitian-definite n by n...
void PTTRS(const OrdinalType n, const OrdinalType nrhs, const ScalarType *d, const ScalarType *e, ScalarType *B, const OrdinalType ldb, OrdinalType *info) const
Solves a tridiagonal system A*X=B using the *D*L' factorization of A computed by PTTRF.
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
Improves the computed solution to a system of linear equations when the coefficient matrix is symmetr...
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
Computes row and column scalings intended to equilibrate an m by n banded matrix A and reduce its con...
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
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
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
Solves an over/underdetermined real m by n linear system A using QR or LQ factorization of A...
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
void LARTG(const ScalarType f, const ScalarType g, MagnitudeType *c, ScalarType *s, ScalarType *r) const
Gnerates a plane rotation that zeros out the second component of the input vector.
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
Solves a system of linear equations A*X=B, where A is a symmetric positive definite matrix factored b...
ScalarType LAPY2(const ScalarType x, const ScalarType y) const
Computes x^2 + y^2 safely, to avoid overflow.
ScalarType LARND(const OrdinalType &idist, OrdinalType *seed) const
void TRTRI(const char UPLO, const char DIAG, const OrdinalType n, const ScalarType *A, const OrdinalType lda, OrdinalType *info) const
Computes the inverse of an upper or lower triangular matrix A.
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
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A...
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
Improves the computed solution to a banded system of linear equations and provides error bounds and b...
void GEBAL(const char JOBZ, const OrdinalType n, ScalarType *A, const OrdinalType lda, OrdinalType ilo, OrdinalType ihi, MagnitudeType *scale, OrdinalType *info) const
Balances a general matrix A, through similarity transformations to make the rows and columns as close...
ScalarTraits< ScalarType >::magnitudeType LANGB(const char NORM, const OrdinalType n, const OrdinalType kl, const OrdinalType ku, const ScalarType *A, const OrdinalType lda, MagnitudeType *WORK) const
Returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of lar...
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
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
Solves a system of linear equations A*X=B or A'*X=B with a general banded n by n matrix A using the L...
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
void POTRF(const char UPLO, const OrdinalType n, ScalarType *A, const OrdinalType lda, OrdinalType *info) const
Computes Cholesky factorization of a real symmetric positive definite matrix A.
Specialization for Sacado::UQ::PCE< Storage<...> >
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
Chooses problem-dependent parameters for the local environment.
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
Solves a system of linear equations A*X=B or A'*X=B or A^H*X=B with a tridiagonal matrix A using the ...
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
Computes the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A in packed stora...
void GESV(const OrdinalType n, const OrdinalType nrhs, ScalarType *A, const OrdinalType lda, OrdinalType *IPIV, ScalarType *B, const OrdinalType ldb, OrdinalType *info) const
Computes the solution to a real system of linear equations A*X=B, where A is factored through GETRF a...
void LATRS(const char UPLO, const char TRANS, const char DIAG, const char NORMIN, const OrdinalType N, ScalarType *A, const OrdinalType LDA, ScalarType *X, MagnitudeType *SCALE, MagnitudeType *CNORM, OrdinalType *INFO) const
Robustly solve a possibly singular triangular linear system.
void GETRF(const OrdinalType m, const OrdinalType n, ScalarType *A, const OrdinalType lda, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a general m by n matrix A using partial pivoting with row interchange...
void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType *A, const OrdinalType lda, ScalarType *B, const OrdinalType ldb, OrdinalType *info) const
Computes the solution to a real system of linear equations A*X=B, where A is a symmetric positive def...
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
Estimates the reciprocal of the condition number of a general real matrix A, in either the 1-norm or ...
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
Reduces a real general matrix A to upper Hessenberg form by orthogonal similarity transformations...
LAPACK(void)
Default Constructor.
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
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
LAPACK(const LAPACK< OrdinalType, ScalarType > &lapack)
Copy Constructor.
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
Estimates the reciprocal of the condition number of a general banded real matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by GETRF.
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
Compute explicit Q factor from QR factorization (GEQRF) (real case).
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, MagnitudeType *FERR, MagnitudeType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Improves the computed solution to a system of linear equations and provides error bounds and backward...
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
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
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
Solves a system of linear equations A*X=B or A'*X=B with a general n by n matrix A using the LU facto...
ScalarType LAMCH(const char CMACH) const
Determines machine parameters for floating point characteristics.
void GEQRF(const OrdinalType m, const OrdinalType n, ScalarType *A, const OrdinalType lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType lwork, OrdinalType *info) const
Computes a QR factorization of a general m by n matrix A.
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
Compute explicit QR factor from QR factorization (GEQRF) (complex case).
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
void GEEQU(const OrdinalType m, const OrdinalType n, const ScalarType *A, const OrdinalType lda, MagnitudeType *R, MagnitudeType *C, MagnitudeType *rowcond, MagnitudeType *colcond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate an m by n matrix A and reduce its condition ...
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
Multiplies the m by n matrix A by the real scalar cto/cfrom.
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
Solves a triangular linear system of the form A*X=B or A**T*X=B, where A is a triangular matrix...
void POEQU(const OrdinalType n, const ScalarType *A, const OrdinalType lda, MagnitudeType *S, MagnitudeType *scond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate a symmetric positive definite matrix A and r...
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
void GTTRF(const OrdinalType n, ScalarType *dl, ScalarType *d, ScalarType *du, ScalarType *du2, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a n by n tridiagonal matrix A using partial pivoting with row interch...
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
Reduces a real symmetric matrix A to tridiagonal form by orthogonal similarity transformations.
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)
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
Generates a real orthogonal matrix Q which is the product of ihi-ilo elementary reflectors of order n...
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
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
Computes for an n by n real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
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
Uses the LU factorization to compute the solution to a real system of linear equations A*X=B...
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
Solves the linear equality-constrained least squares (LSE) problem where A is an m by n matrix...
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