Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Teuchos_LAPACK_UQ_PCE.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Teuchos: Common Tools Package
4 //
5 // Copyright 2004 NTESS and the Teuchos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef _TEUCHOS_LAPACK_UQ_PCE_HPP_
11 #define _TEUCHOS_LAPACK_UQ_PCE_HPP_
12 
13 #include "Teuchos_LAPACK.hpp"
15 #include "Sacado_UQ_PCE.hpp"
16 
24 namespace Teuchos {
25 
26  template<typename OrdinalType, typename Storage>
27  class LAPACK<OrdinalType, Sacado::UQ::PCE<Storage> >
28  {
29  public:
30 
33 
35 
36 
38  inline LAPACK(void) {}
39 
41  inline LAPACK(const LAPACK<OrdinalType, ScalarType>& lapack) {}
42 
44  inline virtual ~LAPACK(void) {}
46 
48 
49 
51  void PTTRF(const OrdinalType n, ScalarType* d, ScalarType* e, OrdinalType* info) const
52  { throw_error("PTTRF"); }
53 
55  void PTTRS(const OrdinalType n, const OrdinalType nrhs, const ScalarType* d, const ScalarType* e, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const
56  { throw_error("PTTRS"); }
57 
59  void POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const
60  { throw_error("POTRF"); }
61 
63  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
64  { throw_error("POTRS"); }
65 
67  void POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const
68  { throw_error("POTRI"); }
69 
71 
72  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
73  { throw_error("POCON"); }
74 
76  void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const
77  { throw_error("POSV"); }
78 
80  void POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, MagnitudeType* S, MagnitudeType* scond, MagnitudeType* amax, OrdinalType* info) const
81  { throw_error("POEQU"); }
82 
84  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
85  { throw_error("PORFS"); }
86 
88  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
89  { throw_error("POSVX"); }
91 
93 
94 
96  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
97  { throw_error("GELS"); }
98 
132  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
133  { throw_error("GELSS"); }
134 
136  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
137  { throw_error("GELSS"); }
138 
140  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
141  { throw_error("GGLSE"); }
142 
144  void GEQRF( const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const
145  { throw_error("GEQRF"); }
146 
148  void GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
149  { throw_error("GETRF"); }
150 
152  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
153  { throw_error("GETRS"); }
154 
156  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
157  { throw_error("LASCL"); }
158 
160  void
161  GEQP3(const OrdinalType m,
162  const OrdinalType n, ScalarType* A,
163  const OrdinalType lda,
164  OrdinalType *jpvt,
165  ScalarType* TAU,
166  ScalarType* WORK,
167  const OrdinalType lwork,
168  MagnitudeType* RWORK,
169  OrdinalType* info ) const
170  { throw_error("GEQP3"); }
171 
173  void
174  LASWP (const OrdinalType N,
175  ScalarType A[],
176  const OrdinalType LDA,
177  const OrdinalType K1,
178  const OrdinalType K2,
179  const OrdinalType IPIV[],
180  const OrdinalType INCX) const
181  { throw_error("LASWP"); }
182 
184  void GBTRF(const OrdinalType m, const OrdinalType n, const OrdinalType kl, const OrdinalType ku, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
185  { throw_error("GBTRF"); }
186 
188  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
189  { throw_error("GBTRS"); }
190 
192  void GTTRF(const OrdinalType n, ScalarType* dl, ScalarType* d, ScalarType* du, ScalarType* du2, OrdinalType* IPIV, OrdinalType* info) const
193  { throw_error("GTTRF"); }
194 
196  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
197  { throw_error("GTTRS"); }
198 
200  void GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const
201  { throw_error("GETRI"); }
202 
207  void
208  LATRS (const char UPLO,
209  const char TRANS,
210  const char DIAG,
211  const char NORMIN,
212  const OrdinalType N,
213  ScalarType* A,
214  const OrdinalType LDA,
215  ScalarType* X,
216  MagnitudeType* SCALE,
217  MagnitudeType* CNORM,
218  OrdinalType* INFO) const
219  { throw_error("LATRS"); }
220 
222  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
223  { throw_error("GECON"); }
224 
226  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
227  { throw_error("GBCON"); }
228 
230  typename 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
231  { throw_error("LANGB"); return MagnitudeType(0); }
232 
234  void GESV(const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const
235  { throw_error("GESV"); }
236 
238  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
239  { throw_error("GEEQU"); }
240 
242  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
243  { throw_error("GERFS"); }
244 
246  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
247  { throw_error("GBEQU"); }
248 
250  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
251  { throw_error("GBRFS"); }
252 
254  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
255  { throw_error("GESVX"); }
256 
260  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
261  { throw_error("SYTRD"); }
262 
264  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
265  { throw_error("GEHRD"); }
266 
268  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
269  { throw_error("TRTRS"); }
270 
272  void TRTRI(const char UPLO, const char DIAG, const OrdinalType n, const ScalarType* A, const OrdinalType lda, OrdinalType* info) const
273  { throw_error("TRTRI"); }
275 
277 
278 
281  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
282  { throw_error("SPEV"); }
283 
287  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
288  { throw_error("SYEV"); }
289 
293  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
294  { throw_error("SYGV"); }
295 
299  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
300  { throw_error("HEEV"); }
301 
305  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
306  { throw_error("HEGV"); }
307 
309  void STEQR(const char COMPZ, const OrdinalType n, ScalarType* D, ScalarType* E, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, OrdinalType* info) const
310  { throw_error("STEQR"); }
312 
314 
315  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
317  { throw_error("HSEQR"); }
318 
322  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
323  { throw_error("GEES"); }
324 
328  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
329  { throw_error("GEES"); }
330 
334  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
335  { throw_error("GEES"); }
336 
342  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
343  { throw_error("GEEV"); }
344 
349  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
350  { throw_error("GEEVX"); }
351 
356  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
357  { throw_error("GGEVX"); }
358 
362  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
363  { throw_error("GGEV"); }
364 
365 
369  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
370  { throw_error("TRSEN"); }
371 
372 
376  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
377  { throw_error("TGSEN"); }
378 
379 
383  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
384  { throw_error("GGES"); }
385 
387 
388 
390 
391  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
393  { throw_error("GESVD"); }
395 
396 
398 
399 
409  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
410  { throw_error("ORMQR"); }
411 
420  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
421  { throw_error("UNMQR"); }
422 
432  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
433  { throw_error("ORGQR"); }
434 
443  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
444  { throw_error("UNGQR"); }
445 
449  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
450  { throw_error("ORGHR"); }
451 
455  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
456  { throw_error("ORMHR"); }
458 
460 
461 
464  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
465  { throw_error("TREVC"); }
466 
470  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
471  { throw_error("TREVC"); }
472 
476  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
477  { throw_error("TREXC"); }
478 
482  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
483  { throw_error("TGEVC"); }
484 
485 
487 
489 
490 
492  void LARTG( const ScalarType f, const ScalarType g, MagnitudeType* c, ScalarType* s, ScalarType* r ) const
493  { throw_error("LARTG"); }
494 
496  void LARFG( const OrdinalType n, ScalarType* alpha, ScalarType* x, const OrdinalType incx, ScalarType* tau ) const
497  { throw_error("LARFG"); }
498 
500 
502 
503 
505  void GEBAL(const char JOBZ, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType ilo, OrdinalType ihi, MagnitudeType* scale, OrdinalType* info) const
506  { throw_error("GEBAL"); }
507 
509  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
510  { throw_error("GEBAK"); }
511 
513 
515 
516  ScalarType LARND( const OrdinalType idist, OrdinalType* seed ) const
518  { throw_error("LARND"); return ScalarType(0); }
519 
521  void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, ScalarType* v ) const
522  { throw_error("LARNV"); }
524 
526 
527 
530  ScalarType LAMCH(const char CMACH) const
531  { throw_error("LAMCH"); return ScalarType(0); }
532 
537  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
538  { throw_error("ILAENV"); return OrdinalType(0); }
540 
542 
543 
546  ScalarType LAPY2(const ScalarType x, const ScalarType y) const
547  { throw_error("LAPY2"); return ScalarType(0); }
549 
550  private:
551 
552  void throw_error(const char *func) const {
554  true, std::logic_error,
555  func << ": Not implemented for Sacado::UQ::PCE scalar type!");
556  }
557 
558  };
559 
560 } // namespace Teuchos
561 
562 #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&#39; 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)
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...
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&#39; 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&#39;*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&lt; Storage&lt;...&gt; &gt;
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&#39;*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...
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&#39;*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