|
Teuchos Package Browser (Single Doxygen Collection)
Version of the Day
|
00001 // @HEADER 00002 // *********************************************************************** 00003 // 00004 // Teuchos: Common Tools Package 00005 // Copyright (2004) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // *********************************************************************** 00040 // @HEADER 00041 00042 #ifndef _TEUCHOS_LAPACK_HPP_ 00043 #define _TEUCHOS_LAPACK_HPP_ 00044 00052 #include "Teuchos_ConfigDefs.hpp" 00053 #include "Teuchos_ScalarTraits.hpp" 00054 00085 namespace Teuchos 00086 { 00087 00088 template<class T> 00089 struct UndefinedLAPACKRoutine 00090 { 00091 // This function should not compile if there is an attempt to instantiate! 00092 static inline T notDefined() { return T::LAPACK_routine_not_defined_for_this_type(); } 00093 }; 00094 00095 template<typename OrdinalType, typename ScalarType> 00096 class LAPACK 00097 { 00098 public: 00099 00100 typedef typename ::Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType; 00101 00103 00104 00106 inline LAPACK(void) {} 00107 00109 inline LAPACK(const LAPACK<OrdinalType, ScalarType>& lapack) {} 00110 00112 inline virtual ~LAPACK(void) {} 00114 00116 00117 00119 void PTTRF(const OrdinalType n, ScalarType* d, ScalarType* e, OrdinalType* info) const; 00120 00122 void PTTRS(const OrdinalType n, const OrdinalType nrhs, const ScalarType* d, const ScalarType* e, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const; 00123 00125 void POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const; 00126 00128 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; 00129 00131 void POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const; 00132 00134 00135 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; 00136 00138 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const; 00139 00141 void POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, MagnitudeType* S, MagnitudeType* scond, MagnitudeType* amax, OrdinalType* info) const; 00142 00144 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; 00145 00147 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; 00149 00151 00152 00154 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; 00155 00189 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; 00190 00192 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; 00193 00195 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; 00196 00198 void GEQRF( const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const; 00199 00201 void GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const; 00202 00204 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; 00205 00207 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; 00208 00210 void 00211 GEQP3(const OrdinalType m, 00212 const OrdinalType n, ScalarType* A, 00213 const OrdinalType lda, 00214 OrdinalType *jpvt, 00215 ScalarType* TAU, 00216 ScalarType* WORK, 00217 const OrdinalType lwork, 00218 MagnitudeType* RWORK, 00219 OrdinalType* info ) const; 00220 00222 void 00223 LASWP (const OrdinalType N, 00224 ScalarType A[], 00225 const OrdinalType LDA, 00226 const OrdinalType K1, 00227 const OrdinalType K2, 00228 const OrdinalType IPIV[], 00229 const OrdinalType INCX) const; 00230 00232 void GBTRF(const OrdinalType m, const OrdinalType n, const OrdinalType kl, const OrdinalType ku, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const; 00233 00235 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; 00236 00238 void GTTRF(const OrdinalType n, ScalarType* dl, ScalarType* d, ScalarType* du, ScalarType* du2, OrdinalType* IPIV, OrdinalType* info) const; 00239 00241 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; 00242 00244 void GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const; 00245 00250 void 00251 LATRS (const char UPLO, 00252 const char TRANS, 00253 const char DIAG, 00254 const char NORMIN, 00255 const OrdinalType N, 00256 ScalarType* A, 00257 const OrdinalType LDA, 00258 ScalarType* X, 00259 MagnitudeType* SCALE, 00260 MagnitudeType* CNORM, 00261 OrdinalType* INFO) const; 00262 00264 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; 00265 00267 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; 00268 00270 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; 00271 00273 void GESV(const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const; 00274 00276 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; 00277 00279 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; 00280 00282 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; 00283 00285 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; 00286 00288 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; 00289 00293 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; 00294 00296 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; 00297 00299 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; 00300 00302 void TRTRI(const char UPLO, const char DIAG, const OrdinalType n, const ScalarType* A, const OrdinalType lda, OrdinalType* info) const; 00304 00306 00307 00310 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; 00311 00315 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; 00316 00320 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; 00321 00325 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; 00326 00330 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; 00331 00333 void STEQR(const char COMPZ, const OrdinalType n, ScalarType* D, ScalarType* E, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, OrdinalType* info) const; 00335 00337 00338 00339 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; 00340 00344 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; 00345 00349 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; 00350 00354 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; 00355 00361 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; 00362 00367 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; 00368 00373 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; 00374 00378 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; 00379 00380 00384 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; 00385 00386 00390 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; 00391 00392 00396 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; 00397 00399 00400 00402 00403 00404 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; 00406 00407 00409 00410 00420 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; 00421 00430 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; 00431 00441 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; 00442 00451 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; 00452 00456 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; 00457 00461 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; 00463 00465 00466 00469 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; 00470 00474 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; 00475 00479 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; 00480 00484 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; 00485 00486 00488 00490 00491 00493 void LARTG( const ScalarType f, const ScalarType g, MagnitudeType* c, ScalarType* s, ScalarType* r ) const; 00494 00496 void LARFG( const OrdinalType n, ScalarType* alpha, ScalarType* x, const OrdinalType incx, ScalarType* tau ) const; 00497 00499 00501 00502 00504 void GEBAL(const char JOBZ, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType ilo, OrdinalType ihi, MagnitudeType* scale, OrdinalType* info) const; 00505 00507 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; 00508 00510 00512 00513 00514 ScalarType LARND( const OrdinalType idist, OrdinalType* seed ) const; 00515 00517 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, ScalarType* v ) const; 00519 00521 00522 00525 ScalarType LAMCH(const char CMACH) const; 00526 00531 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; 00533 00535 00536 00539 ScalarType LAPY2(const ScalarType x, const ScalarType y) const; 00541 }; 00542 00543 // END GENERAL TEMPLATE DECLARATION // 00544 00545 // BEGIN GENERAL TEMPLATE IMPLEMENTATION // 00546 00547 00548 template<typename OrdinalType, typename ScalarType> 00549 void LAPACK<OrdinalType, ScalarType>::PTTRF(const OrdinalType n, ScalarType* d, ScalarType* e, OrdinalType* info) const 00550 { 00551 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00552 } 00553 00554 template<typename OrdinalType, typename ScalarType> 00555 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 00556 { 00557 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00558 } 00559 00560 template<typename OrdinalType, typename ScalarType> 00561 void LAPACK<OrdinalType, ScalarType>::POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const 00562 { 00563 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00564 } 00565 00566 template<typename OrdinalType, typename ScalarType> 00567 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 00568 { 00569 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00570 } 00571 00572 template<typename OrdinalType, typename ScalarType> 00573 void LAPACK<OrdinalType, ScalarType>::POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const 00574 { 00575 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00576 } 00577 00578 template<typename OrdinalType, typename ScalarType> 00579 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 00580 { 00581 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00582 } 00583 00584 template<typename OrdinalType, typename ScalarType> 00585 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 00586 { 00587 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00588 } 00589 00590 template<typename OrdinalType, typename ScalarType> 00591 void LAPACK<OrdinalType, ScalarType>::POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, MagnitudeType* S, MagnitudeType* scond, MagnitudeType* amax, OrdinalType* info) const 00592 { 00593 // Test the input parameters 00594 *info = 0; 00595 if (n < 0) { 00596 *info = -1; 00597 } else if (lda < TEUCHOS_MAX(1, n)) { 00598 *info = -3; 00599 } 00600 if (*info != 0) { 00601 return; 00602 } 00603 00604 ScalarType sZero = ScalarTraits<ScalarType>::zero(); 00605 ScalarType sOne = ScalarTraits<ScalarType>::one(); 00606 MagnitudeType mZero = ScalarTraits<ScalarType>::magnitude(sZero); 00607 MagnitudeType mOne = ScalarTraits<ScalarType>::magnitude(sOne); 00608 00609 // Quick return 00610 if (n == 0) { 00611 *scond = mOne; 00612 *amax = mZero; 00613 return; 00614 } 00615 00616 // Find the minimum and maximum diagonal elements 00617 S[0] = ScalarTraits<ScalarType>::magnitude( A[0] ); 00618 MagnitudeType smin = S[0]; 00619 *amax = S[0]; 00620 for (OrdinalType i=0; i<n; ++i) { 00621 S[i] = ScalarTraits<ScalarType>::magnitude( A[i*lda + i] ); 00622 smin = TEUCHOS_MIN( smin, S[i] ); 00623 *amax = TEUCHOS_MAX( *amax, S[i] ); 00624 } 00625 00626 if (smin < mZero) { 00627 // Find the first non-positve diagonal element and return an error code 00628 for (OrdinalType i=0; i<n; ++i) { 00629 if (S[i] < mZero) 00630 *info = i; 00631 } 00632 } else { 00633 // Set the scale factors to the reciprocals of the diagonal elements 00634 for (OrdinalType i=0; i<n; ++i) { 00635 S[i] = mOne / ScalarTraits<ScalarType>::squareroot( S[i] ); 00636 } 00637 // Compute scond = min(S(i)) / max(S(i)) 00638 *scond = ScalarTraits<ScalarType>::squareroot( smin ) / ScalarTraits<ScalarType>::squareroot( *amax ); 00639 } 00640 } 00641 00642 template<typename OrdinalType, typename ScalarType> 00643 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 00644 { 00645 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00646 } 00647 00648 template<typename OrdinalType, typename ScalarType> 00649 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 00650 { 00651 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00652 } 00653 00654 template<typename OrdinalType, typename ScalarType> 00655 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 00656 { 00657 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00658 } 00659 00660 template<typename OrdinalType, typename ScalarType> 00661 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 00662 { 00663 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00664 } 00665 00666 template<typename OrdinalType, typename ScalarType> 00667 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 00668 { 00669 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00670 } 00671 00672 template<typename OrdinalType, typename ScalarType> 00673 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 00674 { 00675 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00676 } 00677 00678 template<typename OrdinalType, typename ScalarType> 00679 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 00680 { 00681 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00682 } 00683 00684 template<typename OrdinalType, typename ScalarType> 00685 void LAPACK<OrdinalType,ScalarType>::GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const 00686 { 00687 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00688 } 00689 00690 template<typename OrdinalType, typename ScalarType> 00691 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 00692 { 00693 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00694 } 00695 00696 template<typename OrdinalType, typename ScalarType> 00697 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 00698 { 00699 MagnitudeType safeMin = ScalarTraits<ScalarType>::sfmin(); 00700 ScalarType sZero = ScalarTraits<ScalarType>::zero(); 00701 ScalarType sOne = ScalarTraits<ScalarType>::one(); 00702 MagnitudeType mZero = ScalarTraits<ScalarType>::magnitude(sZero); 00703 MagnitudeType mOne = ScalarTraits<ScalarType>::magnitude(sOne); 00704 00705 MagnitudeType smlnum = ScalarTraits<ScalarType>::magnitude(safeMin); 00706 MagnitudeType bignum = ScalarTraits<ScalarType>::magnitude(sOne/smlnum); 00707 00708 OrdinalType i, j; 00709 ScalarType* ptr; 00710 MagnitudeType mul; 00711 bool done = false; 00712 00713 MagnitudeType cfromc = cfrom; 00714 MagnitudeType ctoc = cto; 00715 MagnitudeType cfrom1; 00716 MagnitudeType cto1; 00717 00718 while (!done) { 00719 00720 cfrom1 = cfromc*smlnum; 00721 if (cfrom1 == cfromc) { 00722 // cfromc is an inf. Multiply by a correctly signed zero for finite ctoc, or a NaN if ctoc is infinite. 00723 mul = ctoc / cfromc; 00724 done = true; 00725 cto1 = ctoc; 00726 } else { 00727 cto1 = ctoc / bignum; 00728 if (cto1 == ctoc) { 00729 // ctoc is either 0 or an inf. In both cases, ctoc itself serves as the correct multiplication factor. 00730 mul = ctoc; 00731 done = true; 00732 cfromc = mOne; 00733 } else if (ScalarTraits<ScalarType>::magnitude(cfrom1) > ScalarTraits<ScalarType>::magnitude(ctoc) && ctoc != mZero) { 00734 mul = smlnum; 00735 done = false; 00736 cfromc = cfrom1; 00737 } else if (ScalarTraits<ScalarType>::magnitude(cto1) > ScalarTraits<ScalarType>::magnitude(cfromc)) { 00738 mul = bignum; 00739 done = false; 00740 ctoc = cto1; 00741 } else { 00742 mul = ctoc / cfromc; 00743 done = true; 00744 } 00745 } 00746 00747 for (j=0; j<n; j++) { 00748 ptr = A + j*lda; 00749 for (i=0; i<m; i++) { *ptr = mul * (*ptr); ptr++; } 00750 } 00751 } 00752 00753 } 00754 00755 template<typename OrdinalType, typename ScalarType> 00756 void 00757 LAPACK<OrdinalType,ScalarType>:: 00758 GEQP3 (const OrdinalType m, 00759 const OrdinalType n, 00760 ScalarType* A, 00761 const OrdinalType lda, 00762 OrdinalType *jpvt, 00763 ScalarType* TAU, 00764 ScalarType* WORK, 00765 const OrdinalType lwork, 00766 MagnitudeType* RWORK, 00767 OrdinalType* info) const 00768 { 00769 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00770 } 00771 00772 template<typename OrdinalType, typename ScalarType> 00773 void 00774 LAPACK<OrdinalType, ScalarType>:: 00775 LASWP (const OrdinalType N, 00776 ScalarType A[], 00777 const OrdinalType LDA, 00778 const OrdinalType K1, 00779 const OrdinalType K2, 00780 const OrdinalType IPIV[], 00781 const OrdinalType INCX) const 00782 { 00783 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00784 } 00785 00786 template<typename OrdinalType, typename ScalarType> 00787 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 00788 { 00789 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00790 } 00791 00792 template<typename OrdinalType, typename ScalarType> 00793 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 00794 { 00795 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00796 } 00797 00798 template<typename OrdinalType, typename ScalarType> 00799 void LAPACK<OrdinalType,ScalarType>::GTTRF(const OrdinalType n, ScalarType* dl, ScalarType* d, ScalarType* du, ScalarType* du2, OrdinalType* IPIV, OrdinalType* info) const 00800 { 00801 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00802 } 00803 00804 template<typename OrdinalType, typename ScalarType> 00805 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 00806 { 00807 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00808 } 00809 00810 template<typename OrdinalType, typename ScalarType> 00811 void LAPACK<OrdinalType,ScalarType>::GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const 00812 { 00813 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00814 } 00815 00816 template<typename OrdinalType, typename ScalarType> 00817 void 00818 LAPACK<OrdinalType,ScalarType>:: 00819 LATRS (const char UPLO, 00820 const char TRANS, 00821 const char DIAG, 00822 const char NORMIN, 00823 const OrdinalType N, 00824 ScalarType* A, 00825 const OrdinalType LDA, 00826 ScalarType* X, 00827 MagnitudeType* SCALE, 00828 MagnitudeType* CNORM, 00829 OrdinalType* INFO) const 00830 { 00831 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00832 } 00833 00834 template<typename OrdinalType, typename ScalarType> 00835 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 00836 { 00837 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00838 } 00839 00840 template<typename OrdinalType, typename ScalarType> 00841 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 00842 { 00843 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00844 } 00845 00846 template<typename OrdinalType, typename ScalarType> 00847 typename ScalarTraits<ScalarType>::magnitudeType LAPACK<OrdinalType,ScalarType>::LANGB(const char NORM, const OrdinalType n, const OrdinalType kl, const OrdinalType ku, const ScalarType* A, const OrdinalType lda, MagnitudeType* WORK) const 00848 { 00849 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00850 } 00851 00852 template<typename OrdinalType, typename ScalarType> 00853 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 00854 { 00855 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00856 } 00857 00858 template<typename OrdinalType, typename ScalarType> 00859 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 00860 { 00861 00862 // Test the input parameters 00863 *info = 0; 00864 if (m < 0) { 00865 *info = -1; 00866 } else if (n < 0) { 00867 *info = -2; 00868 } else if (lda < TEUCHOS_MAX(1, m)) { 00869 *info = -4; 00870 } 00871 if (*info != 0) { 00872 return; 00873 } 00874 00875 ScalarType sZero = ScalarTraits<ScalarType>::zero(); 00876 ScalarType sOne = ScalarTraits<ScalarType>::one(); 00877 MagnitudeType mZero = ScalarTraits<ScalarType>::magnitude(sZero); 00878 MagnitudeType mOne = ScalarTraits<ScalarType>::magnitude(sOne); 00879 00880 // Quick return 00881 if (m == 0 || n == 0) { 00882 *rowcond = mOne; 00883 *colcond = mOne; 00884 *amax = mZero; 00885 return; 00886 } 00887 00888 MagnitudeType safeMin = ScalarTraits<ScalarType>::sfmin(); 00889 MagnitudeType smlnum = ScalarTraits<ScalarType>::magnitude(safeMin); 00890 MagnitudeType bignum = ScalarTraits<ScalarType>::magnitude(sOne/smlnum); 00891 00892 // Compute the row scale factors 00893 for (OrdinalType i=0; i<m; i++) { 00894 R[i] = mZero; 00895 } 00896 00897 // Find the maximum element in each row 00898 for (OrdinalType j=0; j<n; j++) { 00899 for (OrdinalType i=0; i<m; i++) { 00900 R[i] = TEUCHOS_MAX( R[i], ScalarTraits<ScalarType>::magnitude( A[j*lda + i] ) ); 00901 } 00902 } 00903 00904 // Find the maximum and minimum scale factors 00905 MagnitudeType rcmin = bignum; 00906 MagnitudeType rcmax = mZero; 00907 for (OrdinalType i=0; i<m; i++) { 00908 rcmax = TEUCHOS_MAX( rcmax, R[i] ); 00909 rcmin = TEUCHOS_MIN( rcmin, R[i] ); 00910 } 00911 *amax = rcmax; 00912 00913 if (rcmin == mZero) { 00914 // Find the first zero scale factor and return an error code 00915 for (OrdinalType i=0; i<m; i++) { 00916 if (R[i] == mZero) 00917 *info = i; 00918 } 00919 } else { 00920 // Invert the scale factors 00921 for (OrdinalType i=0; i<m; i++) { 00922 R[i] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( R[i], smlnum ), bignum ); 00923 } 00924 // Compute rowcond = min(R(i)) / max(R(i)) 00925 *rowcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum ); 00926 } 00927 00928 // Compute the column scale factors 00929 for (OrdinalType j=0; j<n; j++) { 00930 C[j] = mZero; 00931 } 00932 00933 // Find the maximum element in each column, assuming the row scaling computed above 00934 for (OrdinalType j=0; j<n; j++) { 00935 for (OrdinalType i=0; i<m; i++) { 00936 C[j] = TEUCHOS_MAX( C[j], R[i]*ScalarTraits<ScalarType>::magnitude( A[j*lda + i] ) ); 00937 } 00938 } 00939 00940 // Find the maximum and minimum scale factors 00941 rcmin = bignum; 00942 rcmax = mZero; 00943 for (OrdinalType j=0; j<n; j++) { 00944 rcmax = TEUCHOS_MAX( rcmax, C[j] ); 00945 rcmin = TEUCHOS_MIN( rcmin, C[j] ); 00946 } 00947 00948 if (rcmin == mZero) { 00949 // Find the first zero scale factor and return an error code 00950 for (OrdinalType j=0; j<n; j++) { 00951 if (C[j] == mZero) 00952 *info = m+j; 00953 } 00954 } else { 00955 // Invert the scale factors 00956 for (OrdinalType j=0; j<n; j++) { 00957 C[j] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( C[j], smlnum ), bignum ); 00958 } 00959 // Compute colcond = min(C(j)) / max(C(j)) 00960 *colcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum ); 00961 } 00962 } 00963 00964 template<typename OrdinalType, typename ScalarType> 00965 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 00966 { 00967 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 00968 } 00969 00970 template<typename OrdinalType, typename ScalarType> 00971 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 00972 { 00973 00974 // Test the input parameters 00975 * info = 0; 00976 if (m < 0) { 00977 *info = -1; 00978 } else if (n < 0) { 00979 *info = -2; 00980 } else if (kl < 0) { 00981 *info = -3; 00982 } else if (ku < 0) { 00983 *info = -4; 00984 } else if (lda < kl+ku+1) { 00985 *info = -6; 00986 } 00987 if (*info != 0) { 00988 return; 00989 } 00990 00991 ScalarType sZero = ScalarTraits<ScalarType>::zero(); 00992 ScalarType sOne = ScalarTraits<ScalarType>::one(); 00993 MagnitudeType mZero = ScalarTraits<ScalarType>::magnitude(sZero); 00994 MagnitudeType mOne = ScalarTraits<ScalarType>::magnitude(sOne); 00995 00996 // Quick return 00997 if (m == 0 || n == 0) { 00998 *rowcond = mOne; 00999 *colcond = mOne; 01000 *amax = mZero; 01001 return; 01002 } 01003 01004 MagnitudeType safeMin = ScalarTraits<ScalarType>::sfmin(); 01005 MagnitudeType smlnum = ScalarTraits<ScalarType>::magnitude(safeMin); 01006 MagnitudeType bignum = ScalarTraits<ScalarType>::magnitude(sOne/smlnum); 01007 01008 // Compute the row scale factors 01009 for (OrdinalType i=0; i<m; i++) { 01010 R[i] = mZero; 01011 } 01012 01013 // Find the maximum element in each row 01014 for (OrdinalType j=0; j<n; j++) { 01015 for (OrdinalType i=TEUCHOS_MAX(j-ku,0); i<TEUCHOS_MIN(j+kl,m-1); i++) { 01016 R[i] = TEUCHOS_MAX( R[i], ScalarTraits<ScalarType>::magnitude( A[j*lda + ku+i-j] ) ); 01017 } 01018 } 01019 01020 // Find the maximum and minimum scale factors 01021 MagnitudeType rcmin = bignum; 01022 MagnitudeType rcmax = mZero; 01023 for (OrdinalType i=0; i<m; i++) { 01024 rcmax = TEUCHOS_MAX( rcmax, R[i] ); 01025 rcmin = TEUCHOS_MIN( rcmin, R[i] ); 01026 } 01027 *amax = rcmax; 01028 01029 if (rcmin == mZero) { 01030 // Find the first zero scale factor and return an error code 01031 for (OrdinalType i=0; i<m; i++) { 01032 if (R[i] == mZero) 01033 *info = i; 01034 } 01035 } else { 01036 // Invert the scale factors 01037 for (OrdinalType i=0; i<m; i++) { 01038 R[i] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( R[i], smlnum ), bignum ); 01039 } 01040 // Compute rowcond = min(R(i)) / max(R(i)) 01041 *rowcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum ); 01042 } 01043 01044 // Compute the column scale factors 01045 for (OrdinalType j=0; j<n; j++) { 01046 C[j] = mZero; 01047 } 01048 01049 // Find the maximum element in each column, assuming the row scaling computed above 01050 for (OrdinalType j=0; j<n; j++) { 01051 for (OrdinalType i=TEUCHOS_MAX(j-ku,0); i<TEUCHOS_MIN(j+kl,m-1); i++) { 01052 C[j] = TEUCHOS_MAX( C[j], R[i]*ScalarTraits<ScalarType>::magnitude( A[j*lda + ku+i-j] ) ); 01053 } 01054 } 01055 01056 // Find the maximum and minimum scale factors 01057 rcmin = bignum; 01058 rcmax = mZero; 01059 for (OrdinalType j=0; j<n; j++) { 01060 rcmax = TEUCHOS_MAX( rcmax, C[j] ); 01061 rcmin = TEUCHOS_MIN( rcmin, C[j] ); 01062 } 01063 01064 if (rcmin == mZero) { 01065 // Find the first zero scale factor and return an error code 01066 for (OrdinalType j=0; j<n; j++) { 01067 if (C[j] == mZero) 01068 *info = m+j; 01069 } 01070 } else { 01071 // Invert the scale factors 01072 for (OrdinalType j=0; j<n; j++) { 01073 C[j] = mOne / TEUCHOS_MIN( TEUCHOS_MAX( C[j], smlnum ), bignum ); 01074 } 01075 // Compute colcond = min(C(j)) / max(C(j)) 01076 *colcond = TEUCHOS_MAX( rcmin, smlnum ) / TEUCHOS_MIN( rcmax, bignum ); 01077 } 01078 } 01079 01080 template<typename OrdinalType, typename ScalarType> 01081 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 01082 { 01083 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01084 } 01085 01086 template<typename OrdinalType, typename ScalarType> 01087 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 01088 { 01089 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01090 } 01091 01092 template<typename OrdinalType, typename ScalarType> 01093 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 01094 { 01095 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01096 } 01097 01098 template<typename OrdinalType, typename ScalarType> 01099 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 01100 { 01101 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01102 } 01103 01104 template<typename OrdinalType, typename ScalarType> 01105 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 01106 { 01107 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01108 } 01109 01110 template<typename OrdinalType, typename ScalarType> 01111 void LAPACK<OrdinalType,ScalarType>::TRTRI(const char UPLO, const char DIAG, const OrdinalType n, const ScalarType* A, const OrdinalType lda, OrdinalType* info) const 01112 { 01113 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01114 } 01115 01116 template<typename OrdinalType, typename ScalarType> 01117 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 01118 { 01119 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01120 } 01121 01122 template<typename OrdinalType, typename ScalarType> 01123 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 01124 { 01125 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01126 } 01127 01128 template<typename OrdinalType, typename ScalarType> 01129 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 01130 { 01131 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01132 } 01133 01134 template<typename OrdinalType, typename ScalarType> 01135 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 01136 { 01137 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01138 } 01139 01140 template<typename OrdinalType, typename ScalarType> 01141 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 01142 { 01143 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01144 } 01145 01146 template<typename OrdinalType, typename ScalarType> 01147 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 01148 { 01149 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01150 } 01151 01152 template<typename OrdinalType, typename ScalarType> 01153 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 01154 { 01155 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01156 } 01157 01158 template<typename OrdinalType, typename ScalarType> 01159 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 01160 { 01161 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01162 } 01163 01164 template<typename OrdinalType, typename ScalarType> 01165 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 01166 { 01167 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01168 } 01169 01170 template<typename OrdinalType, typename ScalarType> 01171 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 01172 { 01173 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01174 } 01175 01176 template<typename OrdinalType, typename ScalarType> 01177 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 01178 { 01179 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01180 } 01181 01182 template<typename OrdinalType, typename ScalarType> 01183 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 01184 { 01185 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01186 } 01187 01188 template<typename OrdinalType, typename ScalarType> 01189 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 01190 { 01191 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01192 } 01193 01194 template<typename OrdinalType, typename ScalarType> 01195 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 01196 { 01197 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01198 } 01199 01200 template<typename OrdinalType, typename ScalarType> 01201 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 01202 { 01203 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01204 } 01205 01206 01207 template<typename OrdinalType, typename ScalarType> 01208 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 01209 { 01210 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01211 } 01212 01213 01214 template<typename OrdinalType, typename ScalarType> 01215 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 01216 { 01217 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01218 } 01219 01220 01221 template<typename OrdinalType, typename ScalarType> 01222 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 01223 { 01224 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01225 } 01226 01227 template<typename OrdinalType, typename ScalarType> 01228 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 01229 { 01230 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01231 } 01232 01233 template<typename OrdinalType, typename ScalarType> 01234 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 01235 { 01236 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01237 } 01238 01239 template<typename OrdinalType, typename ScalarType> 01240 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 01241 { 01242 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01243 } 01244 01245 template<typename OrdinalType, typename ScalarType> 01246 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 01247 { 01248 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01249 } 01250 01251 template<typename OrdinalType, typename ScalarType> 01252 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 01253 { 01254 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01255 } 01256 01257 template<typename OrdinalType, typename ScalarType> 01258 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 01259 { 01260 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01261 } 01262 01263 template<typename OrdinalType, typename ScalarType> 01264 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 01265 { 01266 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01267 } 01268 01269 template<typename OrdinalType, typename ScalarType> 01270 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 01271 { 01272 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01273 } 01274 01275 template<typename OrdinalType, typename ScalarType> 01276 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 01277 { 01278 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01279 } 01280 01281 01282 template<typename OrdinalType, typename ScalarType> 01283 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 01284 { 01285 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01286 } 01287 01288 01289 template<typename OrdinalType, typename ScalarType> 01290 ScalarType LAPACK<OrdinalType, ScalarType>::LAMCH(const char CMACH) const 01291 { 01292 return UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01293 } 01294 01295 template<typename OrdinalType, typename ScalarType> 01296 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 01297 { 01298 return UndefinedLAPACKRoutine<OrdinalType>::notDefined(); 01299 } 01300 01301 template<typename OrdinalType, typename ScalarType> 01302 ScalarType LAPACK<OrdinalType, ScalarType>::LAPY2(const ScalarType x, const ScalarType y) const 01303 { 01304 return UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01305 } 01306 01307 template<typename OrdinalType, typename ScalarType> 01308 void LAPACK<OrdinalType, ScalarType>::LARTG( const ScalarType f, const ScalarType g, MagnitudeType* c, ScalarType* s, ScalarType* r ) const 01309 { 01310 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01311 } 01312 01313 template<typename OrdinalType, typename ScalarType> 01314 void LAPACK<OrdinalType, ScalarType>::LARFG( const OrdinalType n, ScalarType* alpha, ScalarType* x, const OrdinalType incx, ScalarType* tau ) const 01315 { 01316 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01317 } 01318 01319 template<typename OrdinalType, typename ScalarType> 01320 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 01321 { 01322 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01323 } 01324 01325 01326 template<typename OrdinalType, typename ScalarType> 01327 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 01328 { 01329 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01330 } 01331 01332 template<typename OrdinalType, typename ScalarType> 01333 ScalarType LAPACK<OrdinalType, ScalarType>::LARND( const OrdinalType idist, OrdinalType* seed ) const 01334 { 01335 return UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01336 } 01337 01338 template<typename OrdinalType, typename ScalarType> 01339 void LAPACK<OrdinalType, ScalarType>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, ScalarType* v ) const 01340 { 01341 UndefinedLAPACKRoutine<ScalarType>::notDefined(); 01342 } 01343 01344 // END GENERAL TEMPLATE IMPLEMENTATION // 01345 01346 #ifndef DOXYGEN_SHOULD_SKIP_THIS 01347 01348 // BEGIN INT, FLOAT SPECIALIZATION DECLARATION // 01349 01350 template<> 01351 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, float> 01352 { 01353 public: 01354 inline LAPACK(void) {} 01355 inline LAPACK(const LAPACK<int, float>& lapack) {} 01356 inline virtual ~LAPACK(void) {} 01357 01358 // Symmetric positive definite linear system routines 01359 void POTRF(const char UPLO, const int n, float* A, const int lda, int * info) const; 01360 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; 01361 void PTTRF(const int n, float* d, float* e, int* info) const; 01362 void PTTRS(const int n, const int nrhs, const float* d, const float* e, float* B, const int ldb, int* info) const; 01363 void POTRI(const char UPLO, const int n, float* A, const int lda, int* info) const; 01364 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; 01365 void POSV(const char UPLO, const int n, const int nrhs, float* A, const int lda, float* B, const int ldb, int* info) const; 01366 void POEQU(const int n, const float* A, const int lda, float* S, float* scond, float* amax, int* info) const; 01367 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; 01368 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; 01369 01370 // General Linear System Routines 01371 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; 01372 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; 01373 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; 01374 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; 01375 void GEQRF( const int m, const int n, float* A, const int lda, float* TAU, float* WORK, const int lwork, int* info) const; 01376 void GETRF(const int m, const int n, float* A, const int lda, int* IPIV, int* info) const; 01377 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; 01378 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; 01379 01380 void 01381 GEQP3 (const int m, 01382 const int n, 01383 float* A, 01384 const int lda, 01385 int *jpvt, 01386 float* TAU, 01387 float* WORK, 01388 const int lwork, 01389 float* RWORK, 01390 int* info) const; 01391 01392 void LASWP (const int N, 01393 float A[], 01394 const int LDA, 01395 const int K1, 01396 const int K2, 01397 const int IPIV[], 01398 const int INCX) const; 01399 01400 void GBTRF(const int m, const int n, const int kl, const int ku, float* A, const int lda, int* IPIV, int* info) const; 01401 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; 01402 void GTTRF(const int n, float* dl, float* d, float* du, float* du2, int* IPIV, int* info) const; 01403 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; 01404 01405 01406 void GETRI(const int n, float* A, const int lda, const int* IPIV, float* WORK, const int lwork, int* info) const; 01407 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; 01408 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; 01409 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; 01410 float LANGB(const char NORM, const int n, const int kl, const int ku, const float* A, const int lda, float* WORK) const; 01411 void GESV(const int n, const int nrhs, float* A, const int lda, int* IPIV, float* B, const int ldb, int* info) const; 01412 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; 01413 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; 01414 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; 01415 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; 01416 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; 01417 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; 01418 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; 01419 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; 01420 void TRTRI(const char UPLO, const char DIAG, const int n, const float* A, const int lda, int* info) const; 01421 01422 // Symmetric eigenvalue routines. 01423 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; 01424 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; 01425 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; 01426 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; 01427 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; 01428 void STEQR(const char COMPZ, const int n, float* D, float* E, float* Z, const int ldz, float* WORK, int* info) const; 01429 01430 // Non-Hermitian eigenvalue routines. 01431 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; 01432 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; 01433 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; 01434 01435 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; 01436 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; 01437 01438 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; 01439 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; 01440 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; 01441 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; 01442 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; 01443 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; 01444 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; 01445 01446 // SVD routine 01447 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; 01448 01449 // Orthogonal matrix routines. 01450 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; 01451 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; 01452 01453 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; 01454 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; 01455 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; 01456 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; 01457 01458 // Triangular matrix routines. 01459 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; 01460 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; 01461 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; 01462 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; 01463 01464 // Rotation/reflection generators 01465 void LARTG( const float f, const float g, float* c, float* s, float* r ) const; 01466 void LARFG( const int n, float* alpha, float* x, const int incx, float* tau ) const; 01467 01468 // Matrix balancing routines. 01469 void GEBAL(const char JOBZ, const int n, float* A, const int lda, int ilo, int ihi, float* scale, int* info) const; 01470 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; 01471 01472 // Random number generators 01473 float LARND( const int idist, int* seed ) const; 01474 void LARNV( const int idist, int* seed, const int n, float* v ) const; 01475 01476 // Machine characteristics. 01477 float LAMCH(const char CMACH) const; 01478 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; 01479 01480 // Miscellaneous routines. 01481 float LAPY2(const float x, const float y) const; 01482 01483 }; 01484 01485 // END INT, FLOAT SPECIALIZATION DECLARATION // 01486 01487 // BEGIN INT, DOUBLE SPECIALIZATION DECLARATION // 01488 01489 template<> 01490 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, double> 01491 { 01492 public: 01493 inline LAPACK(void) {} 01494 inline LAPACK(const LAPACK<int, double>& lapack) {} 01495 inline virtual ~LAPACK(void) {} 01496 01497 // Symmetric positive definite linear system routines 01498 void PTTRF(const int n, double* d, double* e, int* info) const; 01499 void PTTRS(const int n, const int nrhs, const double* d, const double* e, double* B, const int ldb, int* info) const; 01500 void POTRF(const char UPLO, const int n, double* A, const int lda, int* info) const; 01501 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; 01502 void POTRI(const char UPLO, const int n, double* A, const int lda, int* info) const; 01503 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; 01504 void POSV(const char UPLO, const int n, const int nrhs, double* A, const int lda, double* B, const int ldb, int* info) const; 01505 void POEQU(const int n, const double* A, const int lda, double* S, double* scond, double* amax, int* info) const; 01506 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; 01507 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; 01508 01509 // General linear system routines 01510 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; 01511 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; 01512 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; 01513 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; 01514 void GEQRF( const int m, const int n, double* A, const int lda, double* TAU, double* WORK, const int lwork, int* info) const; 01515 void GETRF(const int m, const int n, double* A, const int lda, int* IPIV, int* info) const; 01516 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; 01517 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; 01518 01519 void 01520 GEQP3 (const int m, 01521 const int n, 01522 double* A, 01523 const int lda, 01524 int *jpvt, 01525 double* TAU, 01526 double* WORK, 01527 const int lwork, 01528 double* RWORK, 01529 int* info) const; 01530 01531 void LASWP (const int N, 01532 double A[], 01533 const int LDA, 01534 const int K1, 01535 const int K2, 01536 const int IPIV[], 01537 const int INCX) const; 01538 01539 void GBTRF(const int m, const int n, const int kl, const int ku, double* A, const int lda, int* IPIV, int* info) const; 01540 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; 01541 void GTTRF(const int n, double* dl, double* d, double* du, double* du2, int* IPIV, int* info) const; 01542 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; 01543 void GETRI(const int n, double* A, const int lda, const int* IPIV, double* WORK, const int lwork, int* info) const; 01544 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; 01545 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; 01546 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; 01547 double LANGB(const char NORM, const int n, const int kl, const int ku, const double* A, const int lda, double* WORK) const; 01548 void GESV(const int n, const int nrhs, double* A, const int lda, int* IPIV, double* B, const int ldb, int* info) const; 01549 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; 01550 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; 01551 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; 01552 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; 01553 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; 01554 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; 01555 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; 01556 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; 01557 void TRTRI(const char UPLO, const char DIAG, const int n, const double* A, const int lda, int* info) const; 01558 01559 // Symmetric eigenproblem routines. 01560 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; 01561 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; 01562 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; 01563 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; 01564 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; 01565 void STEQR(const char COMPZ, const int n, double* D, double* E, double* Z, const int ldz, double* WORK, int* info) const; 01566 01567 // Non-Hermitian eigenproblem routines. 01568 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; 01569 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; 01570 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; 01571 01572 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; 01573 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; 01574 01575 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; 01576 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; 01577 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; 01578 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; 01579 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; 01580 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; 01581 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; 01582 01583 01584 // SVD routine 01585 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; 01586 01587 // Orthogonal matrix routines. 01588 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; 01589 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; 01590 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; 01591 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; 01592 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; 01593 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; 01594 01595 // Triangular matrix routines. 01596 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; 01597 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; 01598 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; 01599 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; 01600 01601 // Rotation/reflection generators 01602 void LARTG( const double f, const double g, double* c, double* s, double* r ) const; 01603 void LARFG( const int n, double* alpha, double* x, const int incx, double* tau ) const; 01604 01605 // Matrix balancing routines. 01606 void GEBAL(const char JOBZ, const int n, double* A, const int lda, int ilo, int ihi, double* scale, int* info) const; 01607 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; 01608 01609 // Random number generators 01610 double LARND( const int idist, int* seed ) const; 01611 void LARNV( const int idist, int* seed, const int n, double* v ) const; 01612 01613 // Machine characteristic routines. 01614 double LAMCH(const char CMACH) const; 01615 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; 01616 01617 // Miscellaneous routines. 01618 double LAPY2(const double x, const double y) const; 01619 01620 }; 01621 01622 // END INT, DOUBLE SPECIALIZATION DECLARATION // 01623 01624 #ifdef HAVE_TEUCHOS_COMPLEX 01625 01626 // BEGIN INT, COMPLEX<FLOAT> SPECIALIZATION DECLARATION // 01627 01628 template<> 01629 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, std::complex<float> > 01630 { 01631 public: 01632 inline LAPACK(void) {} 01633 inline LAPACK(const LAPACK<int, std::complex<float> >& lapack) {} 01634 inline virtual ~LAPACK(void) {} 01635 01636 // Symmetric positive definite linear system routines 01637 void PTTRF(const int n, std::complex<float>* d, std::complex<float>* e, int* info) const; 01638 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; 01639 void POTRF(const char UPLO, const int n, std::complex<float>* A, const int lda, int* info) const; 01640 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; 01641 void POTRI(const char UPLO, const int n, std::complex<float>* A, const int lda, int* info) const; 01642 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; 01643 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; 01644 void POEQU(const int n, const std::complex<float>* A, const int lda, float* S, float* scond, float* amax, int* info) const; 01645 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; 01646 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; 01647 01648 // General Linear System Routines 01649 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; 01650 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; 01651 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; 01652 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; 01653 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; 01654 void GETRF(const int m, const int n, std::complex<float>* A, const int lda, int* IPIV, int* info) const; 01655 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; 01656 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; 01657 01658 void 01659 GEQP3 (const int m, 01660 const int n, 01661 std::complex<float>* A, 01662 const int lda, 01663 int *jpvt, 01664 std::complex<float>* TAU, 01665 std::complex<float>* WORK, 01666 const int lwork, 01667 float* RWORK, 01668 int* info) const; 01669 01670 void LASWP (const int N, 01671 std::complex<float> A[], 01672 const int LDA, 01673 const int K1, 01674 const int K2, 01675 const int IPIV[], 01676 const int INCX) const; 01677 01678 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; 01679 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; 01680 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; 01681 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; 01682 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; 01683 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; 01684 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; 01685 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; 01686 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; 01687 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; 01688 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; 01689 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; 01690 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; 01691 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; 01692 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; 01693 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; 01694 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; 01695 void TRTRI(const char UPLO, const char DIAG, const int n, const std::complex<float>* A, const int lda, int* info) const; 01696 01697 // Symmetric eigenvalue routines. 01698 void STEQR(const char COMPZ, const int n, float* D, float* E, std::complex<float>* Z, const int ldz, float* WORK, int* info) const; 01699 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; 01700 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; 01701 01702 // Non-Hermitian eigenvalue routines. 01703 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; 01704 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; 01705 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; 01706 01707 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; 01708 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; 01709 01710 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; 01711 01712 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; 01713 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; 01714 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; 01715 01716 // SVD routine 01717 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; 01718 01719 // Triangular matrix routines. 01720 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; 01721 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; 01722 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; 01723 01724 // Rotation/reflection generators 01725 void LARTG( const std::complex<float> f, const std::complex<float> g, float* c, std::complex<float>* s, std::complex<float>* r ) const; 01726 void LARFG( const int n, std::complex<float>* alpha, std::complex<float>* x, const int incx, std::complex<float>* tau ) const; 01727 01728 // Matrix balancing routines. 01729 void GEBAL(const char JOBZ, const int n, std::complex<float>* A, const int lda, int ilo, int ihi, float* scale, int* info) const; 01730 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; 01731 01732 // Random number generators 01733 std::complex<float> LARND( const int idist, int* seed ) const; 01734 void LARNV( const int idist, int* seed, const int n, std::complex<float>* v ) const; 01735 01736 // Machine characteristics 01737 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; 01738 01739 }; 01740 01741 // END INT, COMPLEX<FLOAT> SPECIALIZATION DECLARATION // 01742 01743 // BEGIN INT, COMPLEX<DOUBLE> SPECIALIZATION DECLARATION // 01744 01745 template<> 01746 class TEUCHOSNUMERICS_LIB_DLL_EXPORT LAPACK<int, std::complex<double> > 01747 { 01748 public: 01749 inline LAPACK(void) {} 01750 inline LAPACK(const LAPACK<int, std::complex<double> >& lapack) {} 01751 inline virtual ~LAPACK(void) {} 01752 01753 // Symmetric positive definite linear system routines 01754 void PTTRF(const int n, std::complex<double>* d, std::complex<double>* e, int* info) const; 01755 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; 01756 void POTRF(const char UPLO, const int n, std::complex<double>* A, const int lda, int* info) const; 01757 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; 01758 void POTRI(const char UPLO, const int n, std::complex<double>* A, const int lda, int* info) const; 01759 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; 01760 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; 01761 void POEQU(const int n, const std::complex<double>* A, const int lda, double* S, double* scond, double* amax, int* info) const; 01762 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; 01763 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; 01764 01765 // General Linear System Routines 01766 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; 01767 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; 01768 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; 01769 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; 01770 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; 01771 void GETRF(const int m, const int n, std::complex<double>* A, const int lda, int* IPIV, int* info) const; 01772 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; 01773 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; 01774 01775 void 01776 GEQP3 (const int m, 01777 const int n, 01778 std::complex<double>* A, 01779 const int lda, 01780 int *jpvt, 01781 std::complex<double>* TAU, 01782 std::complex<double>* WORK, 01783 const int lwork, 01784 double* RWORK, 01785 int* info) const; 01786 01787 void LASWP (const int N, 01788 std::complex<double> A[], 01789 const int LDA, 01790 const int K1, 01791 const int K2, 01792 const int IPIV[], 01793 const int INCX) const; 01794 01795 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; 01796 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; 01797 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; 01798 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; 01799 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; 01800 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; 01801 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; 01802 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; 01803 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; 01804 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; 01805 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; 01806 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; 01807 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; 01808 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; 01809 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; 01810 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; 01811 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; 01812 void TRTRI(const char UPLO, const char DIAG, const int n, const std::complex<double>* A, const int lda, int* info) const; 01813 01814 // Symmetric eigenvalue routines. 01815 void STEQR(const char COMPZ, const int n, double* D, double* E, std::complex<double>* Z, const int ldz, double* WORK, int* info) const; 01816 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; 01817 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; 01818 01819 // Non-hermitian eigenvalue routines. 01820 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; 01821 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; 01822 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; 01823 01824 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; 01825 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; 01826 01827 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; 01828 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; 01829 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; 01830 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; 01831 01832 // SVD routine 01833 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; 01834 01835 // Triangular matrix routines. 01836 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; 01837 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; 01838 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; 01839 01840 // Rotation/reflection generators 01841 void LARTG( const std::complex<double> f, const std::complex<double> g, double* c, std::complex<double>* s, std::complex<double>* r ) const; 01842 void LARFG( const int n, std::complex<double>* alpha, std::complex<double>* x, const int incx, std::complex<double>* tau ) const; 01843 01844 // Matrix balancing routines. 01845 void GEBAL(const char JOBZ, const int n, std::complex<double>* A, const int lda, int ilo, int ihi, double* scale, int* info) const; 01846 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; 01847 01848 // Random number generators 01849 std::complex<double> LARND( const int idist, int* seed ) const; 01850 void LARNV( const int idist, int* seed, const int n, std::complex<double>* v ) const; 01851 01852 // Machine characteristics 01853 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; 01854 01855 }; 01856 01857 // END INT, COMPLEX<DOUBLE> SPECIALIZATION DECLARATION // 01858 01859 #endif // HAVE_TEUCHOS_COMPLEX 01860 01861 #endif // DOXYGEN_SHOULD_SKIP_THIS 01862 01863 } // namespace Teuchos 01864 01865 #endif // _TEUCHOS_LAPACK_HPP_
1.7.6.1