00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041 #ifndef _TEUCHOS_LAPACK_HPP_
00042 #define _TEUCHOS_LAPACK_HPP_
00043
00051
00052
00053
00054
00055 #if defined (INTEL_CXML)
00056 unsigned int one=1;
00057 #endif
00058
00059 #ifdef CHAR_MACRO
00060 #undef CHAR_MACRO
00061 #endif
00062 #if defined (INTEL_CXML)
00063 #define CHAR_MACRO(char_var) &char_var, one
00064 #else
00065 #define CHAR_MACRO(char_var) &char_var
00066 #endif
00067
00068 #include "Teuchos_ConfigDefs.hpp"
00069 #include "Teuchos_LAPACK_wrappers.hpp"
00070
00101 namespace Teuchos
00102 {
00103
00104 template<class T>
00105 struct UndefinedLAPACKRoutine
00106 {
00107
00108 static inline T notDefined() { return T::LAPACK_routine_not_defined_for_this_type(); };
00109 };
00110
00111 template<typename OrdinalType, typename ScalarType>
00112 class LAPACK
00113 {
00114 public:
00116
00118 inline LAPACK(void) {};
00119
00121 inline LAPACK(const LAPACK<OrdinalType, ScalarType>& LAPACK) {};
00122
00124 inline virtual ~LAPACK(void) {};
00126
00128
00130 void POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const;
00131
00133 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;
00134
00136 void POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const;
00137
00139
00140 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;
00141
00143 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00144
00146 void POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, ScalarType* S, ScalarType* scond, ScalarType* amax, OrdinalType* info) const;
00147
00149 void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, 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;
00150
00152 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;
00154
00156
00158 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;
00159
00161 void GEQRF( const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00162
00164 void GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
00165
00167 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;
00168
00170 void GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00171
00173 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;
00174
00176 void GESV(const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00177
00179 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;
00180
00182 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;
00183
00185 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;
00186
00190 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;
00191
00193 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;
00195
00197
00200 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;
00201
00205 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;
00206
00210 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;
00211
00213 void STEQR(const char COMPZ, const OrdinalType n, ScalarType* D, ScalarType* E, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, OrdinalType* info) const;
00215
00217
00218 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;
00219
00221 void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, 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;
00222
00224 void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* WR, ScalarType* WI, ScalarType* VL, const OrdinalType ldvl, ScalarType* VR, const OrdinalType ldvr, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00226
00228
00231 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;
00232
00236 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;
00237
00241 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;
00243
00245
00246 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;
00247
00249 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;
00251
00253
00254 ScalarType LARND( const OrdinalType idist, OrdinalType* seed ) const;
00255
00257 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, ScalarType* v ) const;
00259
00261
00264 ScalarType LAMCH(const char CMACH) const;
00266
00268
00271 ScalarType LAPY2(const ScalarType x, const ScalarType y) const;
00273 };
00274
00275
00276
00277
00278
00279
00280 template<typename OrdinalType, typename ScalarType>
00281 void LAPACK<OrdinalType, ScalarType>::POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const
00282 {
00283 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00284 }
00285
00286 template<typename OrdinalType, typename ScalarType>
00287 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
00288 {
00289 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00290 }
00291
00292 template<typename OrdinalType, typename ScalarType>
00293 void LAPACK<OrdinalType, ScalarType>::POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const
00294 {
00295 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00296 }
00297
00298 template<typename OrdinalType, typename ScalarType>
00299 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
00300 {
00301 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00302 }
00303
00304 template<typename OrdinalType, typename ScalarType>
00305 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
00306 {
00307 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00308 }
00309
00310 template<typename OrdinalType, typename ScalarType>
00311 void LAPACK<OrdinalType, ScalarType>::POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, ScalarType* S, ScalarType* scond, ScalarType* amax, OrdinalType* info) const
00312 {
00313 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00314 }
00315
00316 template<typename OrdinalType, typename ScalarType>
00317 void LAPACK<OrdinalType, ScalarType>::PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, 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
00318 {
00319 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00320 }
00321
00322 template<typename OrdinalType, typename ScalarType>
00323 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
00324 {
00325 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00326 }
00327
00328 template<typename OrdinalType, typename ScalarType>
00329 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
00330 {
00331 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00332 }
00333
00334 template<typename OrdinalType, typename ScalarType>
00335 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
00336 {
00337 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00338 }
00339
00340 template<typename OrdinalType, typename ScalarType>
00341 void LAPACK<OrdinalType,ScalarType>::GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
00342 {
00343 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00344 }
00345
00346 template<typename OrdinalType, typename ScalarType>
00347 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
00348 {
00349 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00350 }
00351
00352 template<typename OrdinalType, typename ScalarType>
00353 void LAPACK<OrdinalType,ScalarType>::GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const
00354 {
00355 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00356 }
00357
00358 template<typename OrdinalType, typename ScalarType>
00359 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
00360 {
00361 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00362 }
00363
00364 template<typename OrdinalType, typename ScalarType>
00365 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
00366 {
00367 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00368 }
00369
00370 template<typename OrdinalType, typename ScalarType>
00371 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
00372 {
00373 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00374 }
00375
00376 template<typename OrdinalType, typename ScalarType>
00377 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
00378 {
00379 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00380 }
00381
00382 template<typename OrdinalType, typename ScalarType>
00383 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
00384 {
00385 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00386 }
00387
00388 template<typename OrdinalType, typename ScalarType>
00389 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
00390 {
00391 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00392 }
00393
00394 template<typename OrdinalType, typename ScalarType>
00395 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
00396 {
00397 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00398 }
00399
00400 template<typename OrdinalType, typename ScalarType>
00401 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
00402 {
00403 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00404 }
00405
00406 template<typename OrdinalType, typename ScalarType>
00407 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
00408 {
00409 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00410 }
00411
00412 template<typename OrdinalType, typename ScalarType>
00413 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
00414 {
00415 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00416 }
00417
00418 template<typename OrdinalType, typename ScalarType>
00419 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
00420 {
00421 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00422 }
00423
00424 template<typename OrdinalType, typename ScalarType>
00425 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
00426 {
00427 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00428 }
00429
00430 template<typename OrdinalType, typename ScalarType>
00431 void LAPACK<OrdinalType, ScalarType>::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, 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
00432 {
00433 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00434 }
00435
00436 template<typename OrdinalType, typename ScalarType>
00437 void LAPACK<OrdinalType, ScalarType>::GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* WR, ScalarType* WI, ScalarType* VL, const OrdinalType ldvl, ScalarType* VR, const OrdinalType ldvr, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const
00438 {
00439 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00440 }
00441
00442 template<typename OrdinalType, typename ScalarType>
00443 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
00444 {
00445 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00446 }
00447
00448 template<typename OrdinalType, typename ScalarType>
00449 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
00450 {
00451 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00452 }
00453
00454 template<typename OrdinalType, typename ScalarType>
00455 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
00456 {
00457 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00458 }
00459
00460 template<typename OrdinalType, typename ScalarType>
00461 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
00462 {
00463 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00464 }
00465
00466 template<typename OrdinalType, typename ScalarType>
00467 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
00468 {
00469 UndefinedLAPACKRoutine<ScalarType>::notDefined();
00470 }
00471
00472 template<typename OrdinalType, typename ScalarType>
00473 ScalarType LAPACK<OrdinalType, ScalarType>::LAMCH(const char CMACH) const
00474 {
00475 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
00476 }
00477
00478 template<typename OrdinalType, typename ScalarType>
00479 ScalarType LAPACK<OrdinalType, ScalarType>::LAPY2(const ScalarType x, const ScalarType y) const
00480 {
00481 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
00482 }
00483
00484 template<typename OrdinalType, typename ScalarType>
00485 ScalarType LAPACK<OrdinalType, ScalarType>::LARND( const OrdinalType idist, OrdinalType* seed ) const
00486 {
00487 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
00488 }
00489
00490 template<typename OrdinalType, typename ScalarType>
00491 void LAPACK<OrdinalType, ScalarType>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, ScalarType* v ) const
00492 {
00493 return UndefinedLAPACKRoutine<ScalarType>::notDefined();
00494 }
00495
00496
00497
00498 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00499
00500
00501
00502 template<typename OrdinalType>
00503 class LAPACK<OrdinalType, float>
00504 {
00505 public:
00506 inline LAPACK(void) {};
00507 inline LAPACK(const LAPACK<OrdinalType, float>& LAPACK) {};
00508 inline virtual ~LAPACK(void) {};
00509
00510
00511 void POTRF(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* info) const;
00512 void POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, float* B, const OrdinalType ldb, OrdinalType* info) const;
00513 void POTRI(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* info) const;
00514 void POCON(const char UPLO, const OrdinalType n, const float* A, const OrdinalType lda, const float anorm, float* rcond, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00515 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, OrdinalType* info) const;
00516 void POEQU(const OrdinalType n, const float* A, const OrdinalType lda, float* S, float* scond, float* amax, OrdinalType* info) const;
00517 void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, const float* AF, const OrdinalType ldaf, const float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00518 void POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* AF, const OrdinalType ldaf, char EQUED, float* S, float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00519
00520
00521 void GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00522 void GEQRF( const OrdinalType m, const OrdinalType n, float* A, const OrdinalType lda, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00523 void GETRF(const OrdinalType m, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
00524 void GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, const OrdinalType* IPIV, float* B, const OrdinalType ldb, OrdinalType* info) const;
00525 void GETRI(const OrdinalType n, float* A, const OrdinalType lda, const OrdinalType* IPIV, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00526 void GECON(const char NORM, const OrdinalType n, const float* A, const OrdinalType lda, const float anorm, float* rcond, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00527 void GESV(const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, OrdinalType* IPIV, float* B, const OrdinalType ldb, OrdinalType* info) const;
00528 void GEEQU(const OrdinalType m, const OrdinalType n, const float* A, const OrdinalType lda, float* R, float* C, float* rowcond, float* colcond, float* amax, OrdinalType* info) const;
00529 void GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, const float* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00530 void GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, float* R, float* C, float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00531 void SYTRD(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* D, float* E, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00532 void GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* A, const OrdinalType lda, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00533
00534
00535 void SPEV(const char JOBZ, const char UPLO, const OrdinalType n, float* AP, float* W, float* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const;
00536 void SYEV(const char JOBZ, const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* W, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00537 void SYGV(const OrdinalType itype, const char JOBZ, const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, float* W, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00538 void STEQR(const char COMPZ, const OrdinalType n, float* D, float* E, float* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const;
00539
00540
00541 void HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* H, const OrdinalType ldh, float* WR, float* WI, float* Z, const OrdinalType ldz, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00542 void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* sdim, float* WR, float* WI, float* VS, const OrdinalType ldvs, float* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const;
00543 void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, float* A, const OrdinalType lda, float* WR, float* WI, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00544
00545
00546 void ORMQR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType k, float* A, const OrdinalType lda, const float* TAU, float* C, const OrdinalType ldc, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00547 void ORGHR(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* A, const OrdinalType lda, const float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00548 void ORMHR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, const float* A, const OrdinalType lda, const float* TAU, float* C, const OrdinalType ldc, float* WORK, const OrdinalType lwork, OrdinalType* info) const;
00549
00550
00551 void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const float* T, const OrdinalType ldt, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, float* WORK, OrdinalType* info) const;
00552 void TREXC(const char COMPQ, const OrdinalType n, float* T, const OrdinalType ldt, float* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, float* WORK, OrdinalType* info) const;
00553
00554
00555 float LARND( const OrdinalType idist, OrdinalType* seed ) const;
00556 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, float* v ) const;
00557
00558
00559 float LAMCH(const char CMACH) const;
00560
00561
00562 float LAPY2(const float x, const float y) const;
00563 };
00564
00565
00566
00567
00568
00569 template<typename OrdinalType>
00570 void LAPACK<OrdinalType, float>::POTRF(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* info) const
00571 {
00572 SPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00573 }
00574
00575 template<typename OrdinalType>
00576 void LAPACK<OrdinalType, float>::POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, float* B, const OrdinalType ldb, OrdinalType* info) const
00577 {
00578 SPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00579 }
00580
00581 template<typename OrdinalType>
00582 void LAPACK<OrdinalType, float>::POTRI(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* info) const
00583 {
00584 SPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00585 }
00586
00587 template<typename OrdinalType>
00588 void LAPACK<OrdinalType, float>::POCON(const char UPLO, const OrdinalType n, const float* A, const OrdinalType lda, const float anorm, float* rcond, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00589 {
00590 SPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00591 }
00592
00593 template<typename OrdinalType>
00594 void LAPACK<OrdinalType, float>::POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, OrdinalType* info) const
00595 {
00596 SPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00597 }
00598
00599 template<typename OrdinalType>
00600 void LAPACK<OrdinalType, float>::POEQU(const OrdinalType n, const float* A, const OrdinalType lda, float* S, float* scond, float* amax, OrdinalType* info) const
00601 {
00602 SPOEQU_F77(&n, A, &lda, S, scond, amax, info);
00603 }
00604
00605 template<typename OrdinalType>
00606 void LAPACK<OrdinalType, float>::PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, const float* AF, const OrdinalType ldaf, const float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00607 {
00608 SPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
00609 }
00610
00611 template<typename OrdinalType>
00612 void LAPACK<OrdinalType, float>::POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* AF, const OrdinalType ldaf, char EQUED, float* S, float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00613 {
00614 SPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHAR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info);
00615 }
00616
00617 template<typename OrdinalType>
00618 void LAPACK<OrdinalType,float>::GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00619 {
00620 SGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
00621 }
00622
00623 template<typename OrdinalType>
00624 void LAPACK<OrdinalType,float>::GEQRF( const OrdinalType m, const OrdinalType n, float* A, const OrdinalType lda, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00625 {
00626 SGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
00627 }
00628
00629 template<typename OrdinalType>
00630 void LAPACK<OrdinalType,float>::GETRF(const OrdinalType m, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
00631 {
00632 SGETRF_F77(&m, &n, A, &lda, IPIV, info);
00633 }
00634
00635 template<typename OrdinalType>
00636 void LAPACK<OrdinalType,float>::GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, const OrdinalType* IPIV, float* B, const OrdinalType ldb, OrdinalType* info) const
00637 {
00638 SGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00639 }
00640
00641 template<typename OrdinalType>
00642 void LAPACK<OrdinalType,float>::GETRI(const OrdinalType n, float* A, const OrdinalType lda, const OrdinalType* IPIV, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00643 {
00644 SGETRI_F77(&n, A, &lda, IPIV, WORK, lwork, info);
00645 }
00646
00647 template<typename OrdinalType>
00648 void LAPACK<OrdinalType,float>::GECON(const char NORM, const OrdinalType n, const float* A, const OrdinalType lda, const float anorm, float* rcond, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00649 {
00650 SGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00651 }
00652
00653 template<typename OrdinalType>
00654 void LAPACK<OrdinalType,float>::GESV(const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, OrdinalType* IPIV, float* B, const OrdinalType ldb, OrdinalType* info) const
00655 {
00656 SGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00657 }
00658
00659 template<typename OrdinalType>
00660 void LAPACK<OrdinalType,float>::GEEQU(const OrdinalType m, const OrdinalType n, const float* A, const OrdinalType lda, float* R, float* C, float* rowcond, float* colcond, float* amax, OrdinalType* info) const
00661 {
00662 SGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
00663 }
00664
00665 template<typename OrdinalType>
00666 void LAPACK<OrdinalType,float>::GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const float* A, const OrdinalType lda, const float* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00667 {
00668 SGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
00669 }
00670
00671 template<typename OrdinalType>
00672 void LAPACK<OrdinalType,float>::GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, float* A, const OrdinalType lda, float* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, float* R, float* C, float* B, const OrdinalType ldb, float* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, float* WORK, OrdinalType* IWORK, OrdinalType* info) const
00673 {
00674 SGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHAR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info);
00675 }
00676
00677 template<typename OrdinalType>
00678 void LAPACK<OrdinalType,float>::SYTRD(const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* D, float* E, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00679 {
00680 SSYTRD_F77(CHAR_MACRO(UPLO), &n, A, &lda, D, E, TAU, WORK, &lwork, info);
00681 }
00682
00683 template<typename OrdinalType>
00684 void LAPACK<OrdinalType,float>::GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* A, const OrdinalType lda, float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00685 {
00686 SGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
00687 }
00688
00689 template<typename OrdinalType>
00690 void LAPACK<OrdinalType,float>::SPEV(const char JOBZ, const char UPLO, const OrdinalType n, float* AP, float* W, float* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const
00691 {
00692 SSPEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, AP, W, Z, &ldz, WORK, info);
00693 }
00694
00695 template<typename OrdinalType>
00696 void LAPACK<OrdinalType,float>::SYEV(const char JOBZ, const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* W, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00697 {
00698 SSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info);
00699 }
00700
00701 template<typename OrdinalType>
00702 void LAPACK<OrdinalType,float>::SYGV(const OrdinalType itype, const char JOBZ, const char UPLO, const OrdinalType n, float* A, const OrdinalType lda, float* B, const OrdinalType ldb, float* W, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00703 {
00704 SSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info);
00705 }
00706
00707 template<typename OrdinalType>
00708 void LAPACK<OrdinalType,float>::STEQR(const char COMPZ, const OrdinalType n, float* D, float* E, float* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const
00709 {
00710 SSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
00711 }
00712
00713 template<typename OrdinalType>
00714 void LAPACK<OrdinalType, float>::HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* H, const OrdinalType ldh, float* WR, float* WI, float* Z, const OrdinalType ldz, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00715 {
00716 SHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, WR, WI, Z, &ldz, WORK, &lwork, info);
00717 }
00718
00719 template<typename OrdinalType>
00720 void LAPACK<OrdinalType, float>::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* sdim, float* WR, float* WI, float* VS, const OrdinalType ldvs, float* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const
00721 {
00722 SGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
00723 }
00724
00725 template<typename OrdinalType>
00726 void LAPACK<OrdinalType, float>::GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, float* A, const OrdinalType lda, float* WR, float* WI, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00727 {
00728 SGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, WORK, &lwork, info);
00729 }
00730
00731 template<typename OrdinalType>
00732 void LAPACK<OrdinalType, float>::ORMQR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType k, float* A, const OrdinalType lda, const float* TAU, float* C, const OrdinalType ldc, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00733 {
00734 SORMQR(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
00735 }
00736
00737 template<typename OrdinalType>
00738 void LAPACK<OrdinalType, float>::ORGHR(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, float* A, const OrdinalType lda, const float* TAU, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00739 {
00740 SORGHR_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
00741 }
00742
00743 template<typename OrdinalType>
00744 void LAPACK<OrdinalType, float>::ORMHR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, const float* A, const OrdinalType lda, const float* TAU, float* C, const OrdinalType ldc, float* WORK, const OrdinalType lwork, OrdinalType* info) const
00745 {
00746 SORMHR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &ilo, &ihi, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
00747 }
00748
00749 template<typename OrdinalType>
00750 void LAPACK<OrdinalType, float>::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const float* T, const OrdinalType ldt, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, float* WORK, OrdinalType* info) const
00751 {
00752 if(HOWMNY=='S')
00753 {
00754 *info = -3;
00755 }
00756 else
00757 {
00758 STREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
00759 }
00760 }
00761
00762 template<typename OrdinalType>
00763 void LAPACK<OrdinalType, float>::TREXC(const char COMPQ, const OrdinalType n, float* T, const OrdinalType ldt, float* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, float* WORK, OrdinalType* info) const
00764 {
00765 STREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, WORK, info);
00766 }
00767
00768 template<typename OrdinalType>
00769 float LAPACK<OrdinalType, float>::LARND( const OrdinalType idist, OrdinalType* seed ) const
00770 {
00771 return(SLARND_F77(&idist, seed));
00772 }
00773
00774 template<typename OrdinalType>
00775 void LAPACK<OrdinalType, float>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, float* v ) const
00776 {
00777 SLARNV_F77(&idist, seed, &n, v);
00778 }
00779
00780 template<typename OrdinalType>
00781 float LAPACK<OrdinalType, float>::LAMCH(const char CMACH) const
00782 {
00783 return(SLAMCH_F77(CHAR_MACRO(CMACH)));
00784 }
00785
00786 template<typename OrdinalType>
00787 float LAPACK<OrdinalType, float>::LAPY2(const float x, const float y) const
00788 {
00789 return SLAPY2_F77(&x, &y);
00790 }
00791
00792
00793
00794
00795
00796 template<typename OrdinalType>
00797 class LAPACK<OrdinalType, double>
00798 {
00799 public:
00800 inline LAPACK(void) {};
00801 inline LAPACK(const LAPACK<OrdinalType, double>& LAPACK) {};
00802 inline virtual ~LAPACK(void) {};
00803
00804
00805 void POTRF(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const;
00806 void POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, double* B, const OrdinalType ldb, OrdinalType* info) const;
00807 void POTRI(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const;
00808 void POCON(const char UPLO, const OrdinalType n, const double* A, const OrdinalType lda, const double anorm, double* rcond, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00809 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, OrdinalType* info) const;
00810 void POEQU(const OrdinalType n, const double* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const;
00811 void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, const double* AF, const OrdinalType ldaf, const double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00812 void POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* AF, const OrdinalType ldaf, char EQUED, double* S, double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00813
00814
00815 void GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00816 void GEQRF( const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00817 void GETRF(const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
00818 void GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, const OrdinalType* IPIV, double* B, const OrdinalType ldb, OrdinalType* info) const;
00819 void GETRI(const OrdinalType n, double* A, const OrdinalType lda, const OrdinalType* IPIV, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00820 void GECON(const char NORM, const OrdinalType n, const double* A, const OrdinalType lda, const double anorm, double* rcond, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00821 void GESV(const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, OrdinalType* IPIV, double* B, const OrdinalType ldb, OrdinalType* info) const;
00822 void GEEQU(const OrdinalType m, const OrdinalType n, const double* A, const OrdinalType lda, double* R, double* C, double* rowcond, double* colcond, double* amax, OrdinalType* info) const;
00823 void GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, const double* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00824 void GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, double* R, double* C, double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00825 void SYTRD(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* D, double* E, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00826 void GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* A, const OrdinalType lda, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00827
00828
00829 void SPEV(const char JOBZ, const char UPLO, const OrdinalType n, double* AP, double* W, double* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const;
00830 void SYEV(const char JOBZ, const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* W, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00831 void SYGV(const OrdinalType itype, const char JOBZ, const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, double* W, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00832 void STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, double* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const;
00833
00834
00835 void HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* H, const OrdinalType ldh, double* WR, double* WI, double* Z, const OrdinalType ldz, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00836 void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* sdim, double* WR, double* WI, double* VS, const OrdinalType ldvs, double* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const;
00837 void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, double* A, const OrdinalType lda, double* WR, double* WI, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00838
00839
00840 void ORMQR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType k, double* A, const OrdinalType lda, const double* TAU, double* C, const OrdinalType ldc, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00841 void ORGHR(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* A, const OrdinalType lda, const double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00842 void ORMHR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, const double* A, const OrdinalType lda, const double* TAU, double* C, const OrdinalType ldc, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00843
00844
00845 void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const double* T, const OrdinalType ldt, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, double* WORK, OrdinalType* info) const;
00846 void TREXC(const char COMPQ, const OrdinalType n, double* T, const OrdinalType ldt, double* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, double* WORK, OrdinalType* info) const;
00847
00848
00849 double LARND( const OrdinalType idist, OrdinalType* seed ) const;
00850 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, double* v ) const;
00851
00852
00853 double LAMCH(const char CMACH) const;
00854
00855
00856 double LAPY2(const double x, const double y) const;
00857 };
00858
00859
00860
00861
00862
00863
00864 template<typename OrdinalType>
00865 void LAPACK<OrdinalType, double>::POTRF(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const
00866 {
00867 DPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00868 }
00869
00870 template<typename OrdinalType>
00871 void LAPACK<OrdinalType, double>::POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, double* B, const OrdinalType ldb, OrdinalType* info) const
00872 {
00873 DPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00874 }
00875
00876 template<typename OrdinalType>
00877 void LAPACK<OrdinalType, double>::POTRI(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const
00878 {
00879 DPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00880 }
00881
00882 template<typename OrdinalType>
00883 void LAPACK<OrdinalType, double>::POCON(const char UPLO, const OrdinalType n, const double* A, const OrdinalType lda, const double anorm, double* rcond, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00884 {
00885 DPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00886 }
00887
00888 template<typename OrdinalType>
00889 void LAPACK<OrdinalType, double>::POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, OrdinalType* info) const
00890 {
00891 DPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00892 }
00893
00894 template<typename OrdinalType>
00895 void LAPACK<OrdinalType, double>::POEQU(const OrdinalType n, const double* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const
00896 {
00897 DPOEQU_F77(&n, A, &lda, S, scond, amax, info);
00898 }
00899
00900 template<typename OrdinalType>
00901 void LAPACK<OrdinalType, double>::PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, const double* AF, const OrdinalType ldaf, const double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00902 {
00903 DPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
00904 }
00905
00906 template<typename OrdinalType>
00907 void LAPACK<OrdinalType, double>::POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* AF, const OrdinalType ldaf, char EQUED, double* S, double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00908 {
00909 DPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHAR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info);
00910 }
00911
00912 template<typename OrdinalType>
00913 void LAPACK<OrdinalType,double>::GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00914 {
00915 DGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, &info);
00916 }
00917
00918 template<typename OrdinalType>
00919 void LAPACK<OrdinalType,double>::GEQRF( const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00920 {
00921 DGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
00922 }
00923
00924 template<typename OrdinalType>
00925 void LAPACK<OrdinalType,double>::GETRF(const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
00926 {
00927 DGETRF_F77(&m, &n, A, &lda, IPIV, info);
00928 }
00929
00930 template<typename OrdinalType>
00931 void LAPACK<OrdinalType,double>::GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, const OrdinalType* IPIV, double* B, const OrdinalType ldb, OrdinalType* info) const
00932 {
00933 DGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00934 }
00935
00936 template<typename OrdinalType>
00937 void LAPACK<OrdinalType,double>::GETRI(const OrdinalType n, double* A, const OrdinalType lda, const OrdinalType* IPIV, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00938 {
00939 DGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info);
00940 }
00941
00942 template<typename OrdinalType>
00943 void LAPACK<OrdinalType,double>::GECON(const char NORM, const OrdinalType n, const double* A, const OrdinalType lda, const double anorm, double* rcond, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00944 {
00945 DGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00946 }
00947
00948 template<typename OrdinalType>
00949 void LAPACK<OrdinalType,double>::GESV(const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, OrdinalType* IPIV, double* B, const OrdinalType ldb, OrdinalType* info) const
00950 {
00951 DGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00952 }
00953
00954 template<typename OrdinalType>
00955 void LAPACK<OrdinalType,double>::GEEQU(const OrdinalType m, const OrdinalType n, const double* A, const OrdinalType lda, double* R, double* C, double* rowcond, double* colcond, double* amax, OrdinalType* info) const
00956 {
00957 DGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
00958 }
00959
00960 template<typename OrdinalType>
00961 void LAPACK<OrdinalType,double>::GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const double* A, const OrdinalType lda, const double* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00962 {
00963 DGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
00964 }
00965
00966 template<typename OrdinalType>
00967 void LAPACK<OrdinalType,double>::GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, double* R, double* C, double* B, const OrdinalType ldb, double* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, double* WORK, OrdinalType* IWORK, OrdinalType* info) const
00968 {
00969 DGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHAR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info);
00970 }
00971
00972 template<typename OrdinalType>
00973 void LAPACK<OrdinalType,double>::SYTRD(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* D, double* E, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00974 {
00975 DSYTRD_F77(CHAR_MACRO(UPLO), &n, A, &lda, D, E, TAU, WORK, &lwork, info);
00976 }
00977
00978 template<typename OrdinalType>
00979 void LAPACK<OrdinalType, double>::GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* A, const OrdinalType lda, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00980 {
00981 DGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
00982 }
00983
00984 template<typename OrdinalType>
00985 void LAPACK<OrdinalType,double>::SPEV(const char JOBZ, const char UPLO, const OrdinalType n, double* AP, double* W, double* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const
00986 {
00987 DSPEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, AP, W, Z, &ldz, WORK, info);
00988 }
00989
00990 template<typename OrdinalType>
00991 void LAPACK<OrdinalType,double>::SYEV(const char JOBZ, const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* W, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00992 {
00993 DSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info);
00994 }
00995
00996 template<typename OrdinalType>
00997 void LAPACK<OrdinalType,double>::SYGV(const OrdinalType itype, const char JOBZ, const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, double* W, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00998 {
00999 DSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info);
01000 }
01001
01002 template<typename OrdinalType>
01003 void LAPACK<OrdinalType,double>::STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, double* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const
01004 {
01005 DSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01006 }
01007
01008 template<typename OrdinalType>
01009 void LAPACK<OrdinalType, double>::HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* H, const OrdinalType ldh, double* WR, double* WI, double* Z, const OrdinalType ldz, double* WORK, const OrdinalType lwork, OrdinalType* info) const
01010 {
01011 DHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, WR, WI, Z, &ldz, WORK, &lwork, info);
01012 }
01013
01014 template<typename OrdinalType>
01015 void LAPACK<OrdinalType, double>::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* sdim, double* WR, double* WI, double* VS, const OrdinalType ldvs, double* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const
01016 {
01017 DGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
01018 }
01019
01020 template<typename OrdinalType>
01021 void LAPACK<OrdinalType, double>::GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, double* A, const OrdinalType lda, double* WR, double* WI, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, double* WORK, const OrdinalType lwork, OrdinalType* info) const
01022 {
01023 DGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, WORK, &lwork, info);
01024 }
01025
01026 template<typename OrdinalType>
01027 void LAPACK<OrdinalType, double>::ORMQR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType k, double* A, const OrdinalType lda, const double* TAU, double* C, const OrdinalType ldc, double* WORK, const OrdinalType lwork, OrdinalType* info) const
01028 {
01029 DORMQR(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
01030 }
01031
01032 template<typename OrdinalType>
01033 void LAPACK<OrdinalType, double>::ORGHR(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, double* A, const OrdinalType lda, const double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const
01034 {
01035 DORGHR_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01036 }
01037
01038 template<typename OrdinalType>
01039 void LAPACK<OrdinalType, double>::ORMHR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, const double* A, const OrdinalType lda, const double* TAU, double* C, const OrdinalType ldc, double* WORK, const OrdinalType lwork, OrdinalType* info) const
01040 {
01041 DORMHR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &ilo, &ihi, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
01042 }
01043
01044 template<typename OrdinalType>
01045 void LAPACK<OrdinalType, double>::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const double* T, const OrdinalType ldt, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, double* WORK, OrdinalType* info) const
01046 {
01047 if(HOWMNY=='S')
01048 {
01049 *info = -3;
01050 }
01051 else
01052 {
01053 DTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
01054 }
01055 }
01056
01057 template<typename OrdinalType>
01058 void LAPACK<OrdinalType, double>::TREXC(const char COMPQ, const OrdinalType n, double* T, const OrdinalType ldt, double* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, double* WORK, OrdinalType* info) const
01059 {
01060 DTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, WORK, info);
01061 }
01062
01063 template<typename OrdinalType>
01064 double LAPACK<OrdinalType, double>::LARND( const OrdinalType idist, OrdinalType* seed ) const
01065 {
01066 return(DLARND_F77(&idist, seed));
01067 }
01068
01069 template<typename OrdinalType>
01070 void LAPACK<OrdinalType, double>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, double* v ) const
01071 {
01072 DLARNV_F77(&idist, seed, &n, v);
01073 }
01074
01075 template<typename OrdinalType>
01076 double LAPACK<OrdinalType, double>::LAMCH(const char CMACH) const
01077 {
01078 return(DLAMCH_F77(CHAR_MACRO(CMACH)));
01079 }
01080
01081 template<typename OrdinalType>
01082 double LAPACK<OrdinalType, double>::LAPY2(const double x, const double y) const
01083 {
01084 return DLAPY2_F77(&x, &y);
01085 }
01086
01087
01088
01089 #ifdef HAVE_TEUCHOS_COMPLEX
01090
01091
01092
01093 template<typename OrdinalType>
01094 class LAPACK<OrdinalType, complex<float> >
01095 {
01096 public:
01097 inline LAPACK(void) {};
01098 inline LAPACK(const LAPACK<OrdinalType, complex<float> >& LAPACK) {};
01099 inline virtual ~LAPACK(void) {};
01100
01101
01102 void POTRF(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const;
01103 void POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const;
01104 void POTRI(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const;
01105 void POCON(const char UPLO, const OrdinalType n, const complex<float>* A, const OrdinalType lda, const float anorm, float* rcond, complex<float>* WORK, float* rwork, OrdinalType* info) const;
01106 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const;
01107 void POEQU(const OrdinalType n, const complex<float>* A, const OrdinalType lda, float* S, float* scond, float* amax, OrdinalType* info) const;
01108 void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, const complex<float>* AF, const OrdinalType ldaf, const complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01109 void POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* AF, const OrdinalType ldaf, char EQUED, float* S, complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01110
01111
01112 void GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01113 void GEQRF( const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, complex<float>* TAU, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01114 void GETRF(const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
01115 void GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const;
01116 void GETRI(const OrdinalType n, complex<float>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01117 void GECON(const char NORM, const OrdinalType n, const complex<float>* A, const OrdinalType lda, const float anorm, float* rcond, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01118 void GESV(const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const;
01119 void GEEQU(const OrdinalType m, const OrdinalType n, const complex<float>* A, const OrdinalType lda, float* R, float* C, float* rowcond, float* colcond, float* amax, OrdinalType* info) const;
01120 void GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, const complex<float>* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01121 void GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, float* R, float* C, complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01122 void GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<float>* A, const OrdinalType lda, complex<float>* TAU, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01123
01124
01125 void STEQR(const char COMPZ, const OrdinalType n, float* D, float* E, complex<float>* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const;
01126
01127
01128 void HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<float>* H, const OrdinalType ldh, complex<float>* W, complex<float>* Z, const OrdinalType ldz, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01129 void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* sdim, complex<float>* W, complex<float>* VS, const OrdinalType ldvs, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* BWORK, OrdinalType* info) const;
01130 void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, complex<float>* A, const OrdinalType lda, complex<float>* W, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* info) const;
01131
01132
01133 void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<float>* T, const OrdinalType ldt, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01134 void TREXC(const char COMPQ, const OrdinalType n, complex<float>* T, const OrdinalType ldt, complex<float>* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, OrdinalType* info) const;
01135
01136
01137 complex<float> LARND( const OrdinalType idist, OrdinalType* seed ) const;
01138 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<float>* v ) const;
01139
01140 };
01141
01142
01143
01144
01145
01146 template<typename OrdinalType>
01147 void LAPACK<OrdinalType, complex<float> >::POTRF(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const
01148 {
01149 CPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01150 }
01151
01152 template<typename OrdinalType>
01153 void LAPACK<OrdinalType, complex<float> >::POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const
01154 {
01155 CPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01156 }
01157
01158 template<typename OrdinalType>
01159 void LAPACK<OrdinalType, complex<float> >::POTRI(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const
01160 {
01161 CPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01162 }
01163
01164 template<typename OrdinalType>
01165 void LAPACK<OrdinalType, complex<float> >::POCON(const char UPLO, const OrdinalType n, const complex<float>* A, const OrdinalType lda, const float anorm, float* rcond, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01166 {
01167 CPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01168 }
01169
01170 template<typename OrdinalType>
01171 void LAPACK<OrdinalType, complex<float> >::POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const
01172 {
01173 CPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01174 }
01175
01176 template<typename OrdinalType>
01177 void LAPACK<OrdinalType, complex<float> >::POEQU(const OrdinalType n, const complex<float>* A, const OrdinalType lda, float* S, float* scond, float* amax, OrdinalType* info) const
01178 {
01179 CPOEQU_F77(&n, A, &lda, S, scond, amax, info);
01180 }
01181
01182 template<typename OrdinalType>
01183 void LAPACK<OrdinalType, complex<float> >::PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, const complex<float>* AF, const OrdinalType ldaf, const complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01184 {
01185 CPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01186 }
01187
01188 template<typename OrdinalType>
01189 void LAPACK<OrdinalType, complex<float> >::POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* AF, const OrdinalType ldaf, char EQUED, float* S, complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01190 {
01191 CPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHAR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
01192 }
01193
01194 template<typename OrdinalType>
01195 void LAPACK<OrdinalType,complex<float> >::GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* B, const OrdinalType ldb, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const
01196 {
01197 CGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
01198 }
01199
01200 template<typename OrdinalType>
01201 void LAPACK<OrdinalType,complex<float> >::GEQRF( const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, complex<float>* TAU, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const
01202 {
01203 CGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
01204 }
01205
01206 template<typename OrdinalType>
01207 void LAPACK<OrdinalType,complex<float> >::GETRF(const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
01208 {
01209 CGETRF_F77(&m, &n, A, &lda, IPIV, info);
01210 }
01211
01212 template<typename OrdinalType>
01213 void LAPACK<OrdinalType,complex<float> >::GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const
01214 {
01215 CGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01216 }
01217
01218 template<typename OrdinalType>
01219 void LAPACK<OrdinalType,complex<float> >::GETRI(const OrdinalType n, complex<float>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const
01220 {
01221 CGETRI_F77(&n, A, &lda, IPIV, WORK, lwork, info);
01222 }
01223
01224 template<typename OrdinalType>
01225 void LAPACK<OrdinalType,complex<float> >::GECON(const char NORM, const OrdinalType n, const complex<float>* A, const OrdinalType lda, const float anorm, float* rcond, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01226 {
01227 CGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01228 }
01229
01230 template<typename OrdinalType>
01231 void LAPACK<OrdinalType,complex<float> >::GESV(const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, complex<float>* B, const OrdinalType ldb, OrdinalType* info) const
01232 {
01233 CGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01234 }
01235
01236 template<typename OrdinalType>
01237 void LAPACK<OrdinalType,complex<float> >::GEEQU(const OrdinalType m, const OrdinalType n, const complex<float>* A, const OrdinalType lda, float* R, float* C, float* rowcond, float* colcond, float* amax, OrdinalType* info) const
01238 {
01239 CGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
01240 }
01241
01242 template<typename OrdinalType>
01243 void LAPACK<OrdinalType,complex<float> >::GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<float>* A, const OrdinalType lda, const complex<float>* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01244 {
01245 CGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01246 }
01247
01248 template<typename OrdinalType>
01249 void LAPACK<OrdinalType,complex<float> >::GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, complex<float>* A, const OrdinalType lda, complex<float>* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, float* R, float* C, complex<float>* B, const OrdinalType ldb, complex<float>* X, const OrdinalType ldx, float* rcond, float* FERR, float* BERR, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01250 {
01251 CGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHAR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
01252 }
01253
01254 template<typename OrdinalType>
01255 void LAPACK<OrdinalType,complex<float> >::GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<float>* A, const OrdinalType lda, complex<float>* TAU, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const
01256 {
01257 CGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01258 }
01259
01260 template<typename OrdinalType>
01261 void LAPACK<OrdinalType,complex<float> >::STEQR(const char COMPZ, const OrdinalType n, float* D, float* E, complex<float>* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const
01262 {
01263 CSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01264 }
01265
01266 template<typename OrdinalType>
01267 void LAPACK<OrdinalType, complex<float> >::HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<float>* H, const OrdinalType ldh, complex<float>* W, complex<float>* Z, const OrdinalType ldz, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const
01268 {
01269 CHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
01270 }
01271
01272 template<typename OrdinalType>
01273 void LAPACK<OrdinalType, complex<float> >::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* sdim, complex<float>* W, complex<float>* VS, const OrdinalType ldvs, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* BWORK, OrdinalType* info) const
01274 {
01275 CGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
01276 }
01277
01278 template<typename OrdinalType>
01279 void LAPACK<OrdinalType, complex<float> >::GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, complex<float>* A, const OrdinalType lda, complex<float>* W, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* info) const
01280 {
01281 CGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
01282 }
01283
01284 template<typename OrdinalType>
01285 void LAPACK<OrdinalType, complex<float> >::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<float>* T, const OrdinalType ldt, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01286 {
01287 if(HOWMNY=='S')
01288 {
01289 *info = -3;
01290 }
01291 else
01292 {
01293 CTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
01294 }
01295 }
01296
01297 template<typename OrdinalType>
01298 void LAPACK<OrdinalType, complex<float> >::TREXC(const char COMPQ, const OrdinalType n, complex<float>* T, const OrdinalType ldt, complex<float>* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, OrdinalType* info) const
01299 {
01300 CTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, info);
01301 }
01302
01303 template<typename OrdinalType>
01304 complex<float> LAPACK<OrdinalType, complex<float> >::LARND( const OrdinalType idist, OrdinalType* seed ) const
01305 {
01306 return(CLARND_F77(&idist, seed));
01307 }
01308
01309 template<typename OrdinalType>
01310 void LAPACK<OrdinalType, complex<float> >::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<float>* v ) const
01311 {
01312 CLARNV_F77(&idist, seed, &n, v);
01313 }
01314
01315
01316
01317
01318
01319 template<typename OrdinalType>
01320 class LAPACK<OrdinalType, complex<double> >
01321 {
01322 public:
01323 inline LAPACK(void) {};
01324 inline LAPACK(const LAPACK<OrdinalType, complex<double> >& LAPACK) {};
01325 inline virtual ~LAPACK(void) {};
01326
01327
01328 void POTRF(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const;
01329 void POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const;
01330 void POTRI(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const;
01331 void POCON(const char UPLO, const OrdinalType n, const complex<double>* A, const OrdinalType lda, const double anorm, double* rcond, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01332 void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const;
01333 void POEQU(const OrdinalType n, const complex<double>* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const;
01334 void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, const complex<double>* AF, const OrdinalType ldaf, const complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01335 void POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* AF, const OrdinalType ldaf, char EQUED, double* S, complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01336
01337
01338 void GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01339 void GEQRF( const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, complex<double>* TAU, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01340 void GETRF(const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
01341 void GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const;
01342 void GETRI(const OrdinalType n, complex<double>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01343 void GECON(const char NORM, const OrdinalType n, const complex<double>* A, const OrdinalType lda, const double anorm, double* rcond, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01344 void GESV(const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const;
01345 void GEEQU(const OrdinalType m, const OrdinalType n, const complex<double>* A, const OrdinalType lda, double* R, double* C, double* rowcond, double* colcond, double* amax, OrdinalType* info) const;
01346 void GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, const complex<double>* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01347 void GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, double* R, double* C, complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01348 void GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<double>* A, const OrdinalType lda, complex<double>* TAU, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01349
01350
01351 void STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, complex<double>* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const;
01352
01353
01354 void HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<double>* H, const OrdinalType ldh, complex<double>* W, complex<double>* Z, const OrdinalType ldz, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01355 void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* sdim, complex<double>* W, complex<double>* VS, const OrdinalType ldvs, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* BWORK, OrdinalType* info) const;
01356 void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, complex<double>* A, const OrdinalType lda, complex<double>* W, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* info) const;
01357
01358
01359 void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<double>* T, const OrdinalType ldt, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01360 void TREXC(const char COMPQ, const OrdinalType n, complex<double>* T, const OrdinalType ldt, complex<double>* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, OrdinalType* info) const;
01361
01362
01363 complex<double> LARND( const OrdinalType idist, OrdinalType* seed ) const;
01364 void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<double>* v ) const;
01365
01366 };
01367
01368
01369
01370
01371
01372 template<typename OrdinalType>
01373 void LAPACK<OrdinalType, complex<double> >::POTRF(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const
01374 {
01375 ZPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01376 }
01377
01378 template<typename OrdinalType>
01379 void LAPACK<OrdinalType, complex<double> >::POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const
01380 {
01381 ZPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01382 }
01383
01384 template<typename OrdinalType>
01385 void LAPACK<OrdinalType, complex<double> >::POTRI(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const
01386 {
01387 ZPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01388 }
01389
01390 template<typename OrdinalType>
01391 void LAPACK<OrdinalType, complex<double> >::POCON(const char UPLO, const OrdinalType n, const complex<double>* A, const OrdinalType lda, const double anorm, double* rcond, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01392 {
01393 ZPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01394 }
01395
01396 template<typename OrdinalType>
01397 void LAPACK<OrdinalType, complex<double> >::POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const
01398 {
01399 ZPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01400 }
01401
01402 template<typename OrdinalType>
01403 void LAPACK<OrdinalType, complex<double> >::POEQU(const OrdinalType n, const complex<double>* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const
01404 {
01405 ZPOEQU_F77(&n, A, &lda, S, scond, amax, info);
01406 }
01407
01408 template<typename OrdinalType>
01409 void LAPACK<OrdinalType, complex<double> >::PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, const complex<double>* AF, const OrdinalType ldaf, const complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01410 {
01411 ZPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01412 }
01413
01414 template<typename OrdinalType>
01415 void LAPACK<OrdinalType, complex<double> >::POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* AF, const OrdinalType ldaf, char EQUED, double* S, complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01416 {
01417 ZPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHAR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
01418 }
01419
01420 template<typename OrdinalType>
01421 void LAPACK<OrdinalType,complex<double> >::GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* B, const OrdinalType ldb, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const
01422 {
01423 ZGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
01424 }
01425
01426 template<typename OrdinalType>
01427 void LAPACK<OrdinalType,complex<double> >::GEQRF( const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, complex<double>* TAU, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const
01428 {
01429 ZGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
01430 }
01431
01432 template<typename OrdinalType>
01433 void LAPACK<OrdinalType,complex<double> >::GETRF(const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
01434 {
01435 ZGETRF_F77(&m, &n, A, &lda, IPIV, info);
01436 }
01437
01438 template<typename OrdinalType>
01439 void LAPACK<OrdinalType,complex<double> >::GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const
01440 {
01441 ZGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01442 }
01443
01444 template<typename OrdinalType>
01445 void LAPACK<OrdinalType,complex<double> >::GETRI(const OrdinalType n, complex<double>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const
01446 {
01447 ZGETRI_F77(&n, A, &lda, IPIV, WORK, lwork, info);
01448 }
01449
01450 template<typename OrdinalType>
01451 void LAPACK<OrdinalType,complex<double> >::GECON(const char NORM, const OrdinalType n, const complex<double>* A, const OrdinalType lda, const double anorm, double* rcond, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01452 {
01453 ZGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01454 }
01455
01456 template<typename OrdinalType>
01457 void LAPACK<OrdinalType,complex<double> >::GESV(const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, complex<double>* B, const OrdinalType ldb, OrdinalType* info) const
01458 {
01459 ZGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01460 }
01461
01462 template<typename OrdinalType>
01463 void LAPACK<OrdinalType,complex<double> >::GEEQU(const OrdinalType m, const OrdinalType n, const complex<double>* A, const OrdinalType lda, double* R, double* C, double* rowcond, double* colcond, double* amax, OrdinalType* info) const
01464 {
01465 ZGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
01466 }
01467
01468 template<typename OrdinalType>
01469 void LAPACK<OrdinalType,complex<double> >::GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const complex<double>* A, const OrdinalType lda, const complex<double>* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01470 {
01471 ZGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01472 }
01473
01474 template<typename OrdinalType>
01475 void LAPACK<OrdinalType,complex<double> >::GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, complex<double>* A, const OrdinalType lda, complex<double>* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, double* R, double* C, complex<double>* B, const OrdinalType ldb, complex<double>* X, const OrdinalType ldx, double* rcond, double* FERR, double* BERR, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01476 {
01477 ZGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHAR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
01478 }
01479
01480 template<typename OrdinalType>
01481 void LAPACK<OrdinalType,complex<double> >::GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<double>* A, const OrdinalType lda, complex<double>* TAU, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const
01482 {
01483 ZGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01484 }
01485
01486 template<typename OrdinalType>
01487 void LAPACK<OrdinalType,complex<double> >::STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, complex<double>* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const
01488 {
01489 ZSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01490 }
01491
01492 template<typename OrdinalType>
01493 void LAPACK<OrdinalType, complex<double> >::HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, complex<double>* H, const OrdinalType ldh, complex<double>* W, complex<double>* Z, const OrdinalType ldz, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const
01494 {
01495 ZHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
01496 }
01497
01498 template<typename OrdinalType>
01499 void LAPACK<OrdinalType, complex<double> >::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* sdim, complex<double>* W, complex<double>* VS, const OrdinalType ldvs, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* BWORK, OrdinalType* info) const
01500 {
01501 ZGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
01502 }
01503
01504 template<typename OrdinalType>
01505 void LAPACK<OrdinalType, complex<double> >::GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, complex<double>* A, const OrdinalType lda, complex<double>* W, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* info) const
01506 {
01507 ZGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
01508 }
01509
01510 template<typename OrdinalType>
01511 void LAPACK<OrdinalType, complex<double> >::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<double>* T, const OrdinalType ldt, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01512 {
01513 if(HOWMNY=='S')
01514 {
01515 *info = -3;
01516 }
01517 else
01518 {
01519 ZTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
01520 }
01521 }
01522
01523 template<typename OrdinalType>
01524 void LAPACK<OrdinalType, complex<double> >::TREXC(const char COMPQ, const OrdinalType n, complex<double>* T, const OrdinalType ldt, complex<double>* Q, const OrdinalType ldq, OrdinalType ifst, OrdinalType ilst, OrdinalType* info) const
01525 {
01526 ZTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, info);
01527 }
01528
01529 template<typename OrdinalType>
01530 complex<double> LAPACK<OrdinalType, complex<double> >::LARND( const OrdinalType idist, OrdinalType* seed ) const
01531 {
01532 return(ZLARND_F77(&idist, seed));
01533 }
01534
01535 template<typename OrdinalType>
01536 void LAPACK<OrdinalType, complex<double> >::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<double>* v ) const
01537 {
01538 ZLARNV_F77(&idist, seed, &n, v);
01539 }
01540
01541
01542
01543 #endif
01544
01545 #endif
01546
01547 }
01548
01549 #endif // _TEUCHOS_LAPACK_HPP_