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