Teuchos_LAPACK.hpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 // 
00004 //                    Teuchos: Common Tools Package
00005 //                 Copyright (2004) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 
00025 // 
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 
00030 // Kris
00031 // 06.11.03 -- Format cleanup
00032 // 06.13.03 -- Initial templatization (Tpetra_LAPACK.cpp is no longer needed)
00033 // 06.17.03 -- Changed LAPACK wrapper calls from XXXXX_F77() to xxxxx_()
00034 //          -- Added warning messages to default function calls
00035 //          -- Added LAPY2 and GEES by request
00036 // 06.18.03 -- Renamed GEES parameters for default implementation
00037 //          -- Changed LAPACK wrapper calls back to XXXXX_F77() from xxxxx_() (Oops!)
00038 //          -- Streamlined character macro stuff
00039 // 07.08.03 -- Move into Teuchos package/namespace
00040 
00041 #ifndef _TEUCHOS_LAPACK_HPP_
00042 #define _TEUCHOS_LAPACK_HPP_
00043 
00051 /* for INTEL_CXML, the second arg may need to be changed to 'one'.  If so
00052 the appropriate declaration of one will need to be added back into
00053 functions that include the macro:
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     // This function should not compile if there is an attempt to instantiate!
00108     static inline T notDefined() { return T::LAPACK_routine_not_defined_for_this_type(); };
00109   };
00110 
00111   template<typename OrdinalType, typename ScalarType>
00112   class LAPACK
00113   {    
00114   public:
00116 
00118     inline LAPACK(void) {};
00119 
00121     inline LAPACK(const LAPACK<OrdinalType, ScalarType>& LAPACK) {};
00122 
00124     inline virtual ~LAPACK(void) {};
00126 
00128 
00130     void POTRF(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const;
00131 
00133     void POTRS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, const ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00134 
00136     void POTRI(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* info) const;
00137 
00139 
00140     void POCON(const char UPLO, const OrdinalType n, const ScalarType* A, const OrdinalType lda, const ScalarType anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00141 
00143     void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00144 
00146     void POEQU(const OrdinalType n, const ScalarType* A, const OrdinalType lda, ScalarType* S, ScalarType* scond, ScalarType* amax, OrdinalType* info) const;
00147 
00149     void PORFS(const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, const ScalarType* AF, const OrdinalType ldaf, const ScalarType* B, const OrdinalType ldb, ScalarType* X, const OrdinalType ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00150 
00152     void POSVX(const char FACT, const char UPLO, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* AF, const OrdinalType ldaf, char EQUED, ScalarType* S, ScalarType* B, const OrdinalType ldb, ScalarType* X, const OrdinalType ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00154 
00156 
00158     void GELS(const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00159 
00161     void GEQRF( const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00162 
00164     void GETRF(const OrdinalType m, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
00165 
00167     void GETRS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00168 
00170     void GETRI(const OrdinalType n, ScalarType* A, const OrdinalType lda, const OrdinalType* IPIV, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00171 
00173     void GECON(const char NORM, const OrdinalType n, const ScalarType* A, const OrdinalType lda, const ScalarType anorm, ScalarType* rcond, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00174 
00176     void GESV(const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, OrdinalType* IPIV, ScalarType* B, const OrdinalType ldb, OrdinalType* info) const;
00177 
00179     void GEEQU(const OrdinalType m, const OrdinalType n, const ScalarType* A, const OrdinalType lda, ScalarType* R, ScalarType* C, ScalarType* rowcond, ScalarType* colcond, ScalarType* amax, OrdinalType* info) const;
00180 
00182     void GERFS(const char TRANS, const OrdinalType n, const OrdinalType nrhs, const ScalarType* A, const OrdinalType lda, const ScalarType* AF, const OrdinalType ldaf, const OrdinalType* IPIV, const ScalarType* B, const OrdinalType ldb, ScalarType* X, const OrdinalType ldx, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00183 
00185     void GESVX(const char FACT, const char TRANS, const OrdinalType n, const OrdinalType nrhs, ScalarType* A, const OrdinalType lda, ScalarType* AF, const OrdinalType ldaf, OrdinalType* IPIV, char EQUED, ScalarType* R, ScalarType* C, ScalarType* B, const OrdinalType ldb, ScalarType* X, const OrdinalType ldx, ScalarType* rcond, ScalarType* FERR, ScalarType* BERR, ScalarType* WORK, OrdinalType* IWORK, OrdinalType* info) const;
00186 
00190     void SYTRD(const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* D, ScalarType* E, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00191 
00193     void GEHRD(const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, ScalarType* A, const OrdinalType lda, ScalarType* TAU, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00195 
00197 
00200     void SPEV(const char JOBZ, const char UPLO, const OrdinalType n, ScalarType* AP, ScalarType* W, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, OrdinalType* info) const;
00201 
00205     void SYEV(const char JOBZ, const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* W, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00206 
00210     void SYGV(const OrdinalType itype, const char JOBZ, const char UPLO, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* B, const OrdinalType ldb, ScalarType* W, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00211 
00213     void STEQR(const char COMPZ, const OrdinalType n, ScalarType* D, ScalarType* E, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, OrdinalType* info) const;
00215 
00217 
00218     void HSEQR(const char JOB, const char COMPZ, const OrdinalType n, const OrdinalType ilo, const OrdinalType ihi, ScalarType* H, const OrdinalType ldh, ScalarType* WR, ScalarType* WI, ScalarType* Z, const OrdinalType ldz, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00219     
00221     void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* sdim, ScalarType* WR, ScalarType* WI, ScalarType* VS, const OrdinalType ldvs, ScalarType* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const;    
00222 
00224     void GEEV(const char JOBVL, const char JOBVR, const OrdinalType n, ScalarType* A, const OrdinalType lda, ScalarType* WR, ScalarType* WI, ScalarType* VL, const OrdinalType ldvl, ScalarType* VR, const OrdinalType ldvr, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00226 
00228 
00231     void ORMQR(const char SIDE, const char TRANS, const OrdinalType m, const OrdinalType n, const OrdinalType k, ScalarType* A, const OrdinalType lda, const ScalarType* TAU, ScalarType* C, const OrdinalType ldc, ScalarType* WORK, const OrdinalType lwork, OrdinalType* info) const;
00232 
00236     void 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* SELECT, const OrdinalType n, const ScalarType* T, const OrdinalType ldt, ScalarType* VL, const OrdinalType ldvl, ScalarType* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, ScalarType* WORK, OrdinalType* info) const;
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     int ILAENV( const int ispec, const string NAME, const string OPTS, const int N1 = -1, const int N2 = -1, const int N3 = -1, const int N4 = -1 ) const;
00277 
00279 
00282     ScalarType LAPY2(const ScalarType x, const ScalarType y) const;
00284   };
00285 
00286   // END GENERAL TEMPLATE DECLARATION //
00287 
00288   // BEGIN GENERAL TEMPLATE IMPLEMENTATION //
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* SELECT, const OrdinalType n, ScalarType* A, const OrdinalType lda, OrdinalType* sdim, ScalarType* WR, ScalarType* WI, ScalarType* VS, const OrdinalType ldvs, ScalarType* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const    
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* SELECT, const OrdinalType n, const ScalarType* T, const OrdinalType ldt, ScalarType* VL, const OrdinalType ldvl, ScalarType* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, ScalarType* WORK, OrdinalType* info) const
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   int LAPACK<OrdinalType, ScalarType>::ILAENV( const int ispec, const string NAME, const string OPTS, const int N1, const int N2, const int N3, const int 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   // END GENERAL TEMPLATE IMPLEMENTATION //
00520 
00521 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00522 
00523   // BEGIN FLOAT PARTIAL SPECIALIZATION DECLARATION //
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     // Symmetric positive definite linear system routines
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     // General Linear System Routines
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     // Symmetric eigenvalue routines.
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     // Hessenberg eigenvalue routines.
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* SELECT, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* sdim, float* WR, float* WI, float* VS, const OrdinalType ldvs, float* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const;    
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     // Orthogonal matrix routines.
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     // Triangular matrix routines.
00575     void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const float* T, const OrdinalType ldt, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, float* WORK, OrdinalType* info) const;
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     // Random number generators
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     // Machine characteristics.
00583     float LAMCH(const char CMACH) const;
00584     int ILAENV( const int ispec, const string NAME, const string OPTS, const int N1 = -1, const int N2 = -1, const int N3 = -1, const int N4 = -1 ) const;
00585 
00586 
00587     // Miscellaneous routines.
00588     float LAPY2(const float x, const float y) const;
00589   };
00590 
00591   // END FLOAT PARTIAL SPECIALIZATION DECLARATION //
00592 
00593   // BEGIN FLOAT PARTIAL SPECIALIZATION IMPLEMENTATION //
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* SELECT, const OrdinalType n, float* A, const OrdinalType lda, OrdinalType* sdim, float* WR, float* WI, float* VS, const OrdinalType ldvs, float* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const    
00747   {
00748     SGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &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* SELECT, const OrdinalType n, const float* T, const OrdinalType ldt, float* VL, const OrdinalType ldvl, float* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, float* WORK, OrdinalType* info) const
00783   {
00784     if(HOWMNY=='S')
00785       {
00786   *info = -3; // We do not support 'S' since it requires a logical function (yuck!)
00787       }
00788     else
00789       {
00790   STREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
00791       }
00792   }
00793   
00794   template<typename OrdinalType>
00795   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
00796   {
00797     STREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, WORK, info);
00798   }
00799 
00800   template<typename OrdinalType>
00801   float LAPACK<OrdinalType, float>::LARND( const OrdinalType idist, OrdinalType* seed ) const
00802   {
00803     return(SLARND_F77(&idist, seed));
00804   }
00805 
00806   template<typename OrdinalType>
00807   void LAPACK<OrdinalType, float>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, float* v ) const
00808   {
00809     SLARNV_F77(&idist, seed, &n, v);
00810   }
00811 
00812   template<typename OrdinalType>
00813   float LAPACK<OrdinalType, float>::LAMCH(const char CMACH) const
00814   {
00815     return(SLAMCH_F77(CHAR_MACRO(CMACH)));
00816   }
00817 
00818   template<typename OrdinalType>
00819   int LAPACK<OrdinalType, float>::ILAENV( const int ispec, const string NAME, const string OPTS, const int N1, const int N2, const int N3, const int N4 ) const
00820   {
00821     unsigned int opts_length = OPTS.length();
00822     string temp_NAME = "s" + NAME;
00823     unsigned int name_length = temp_NAME.length();
00824 #if defined (INTEL_CXML)
00825     return ILAENV_F77(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4 );
00826 #else
00827     return ILAENV_F77(&ispec, &temp_NAME[0], &OPTS[0], &N1, &N2, &N3, &N4, name_length, opts_length );
00828 #endif
00829   }
00830  
00831   template<typename OrdinalType>
00832   float LAPACK<OrdinalType, float>::LAPY2(const float x, const float y) const
00833   {
00834     return SLAPY2_F77(&x, &y);
00835   }
00836 
00837   // END FLOAT PARTIAL SPECIALIZATION IMPLEMENTATION //
00838 
00839   // BEGIN DOUBLE PARTIAL SPECIALIZATION DECLARATION //
00840 
00841   template<typename OrdinalType>
00842   class LAPACK<OrdinalType, double>
00843   {    
00844   public:
00845     inline LAPACK(void) {};
00846     inline LAPACK(const LAPACK<OrdinalType, double>& LAPACK) {};
00847     inline virtual ~LAPACK(void) {};
00848 
00849     // Symmetric positive definite linear system routines
00850     void POTRF(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const;
00851     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;
00852     void POTRI(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const;
00853     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;
00854     void POSV(const char UPLO, const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, double* B, const OrdinalType ldb, OrdinalType* info) const;
00855     void POEQU(const OrdinalType n, const double* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const;
00856     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;
00857     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; 
00858 
00859     // General linear system routines
00860     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;
00861     void GEQRF( const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, double* TAU, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00862     void GETRF(const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
00863     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;
00864     void GETRI(const OrdinalType n, double* A, const OrdinalType lda, const OrdinalType* IPIV, double* WORK, const OrdinalType lwork, OrdinalType* info) const;
00865     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;
00866     void GESV(const OrdinalType n, const OrdinalType nrhs, double* A, const OrdinalType lda, OrdinalType* IPIV, double* B, const OrdinalType ldb, OrdinalType* info) const;
00867     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;
00868     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;
00869     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;
00870     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;
00871     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;
00872 
00873     // Symmetric eigenproblem routines.
00874     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;
00875     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;
00876     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;
00877     void STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, double* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const;
00878 
00879     // Hessenberg eigenproblem routines.
00880     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;
00881     void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* sdim, double* WR, double* WI, double* VS, const OrdinalType ldvs, double* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const;    
00882     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;
00883 
00884     // Orthogonal matrix routines.
00885     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;
00886     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;
00887     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;
00888     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;
00889 
00890     // Triangular matrix routines.
00891     void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const double* T, const OrdinalType ldt, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, double* WORK, OrdinalType* info) const;
00892     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;
00893 
00894     // Random number generators
00895     double LARND( const OrdinalType idist, OrdinalType* seed ) const;
00896     void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, double* v ) const;    
00897 
00898     // Machine characteristic routines.
00899     double LAMCH(const char CMACH) const;
00900     int ILAENV( const int ispec, const string NAME, const string OPTS, const int N1 = -1, const int N2 = -1, const int N3 = -1, const int N4 = -1 ) const;
00901 
00902     // Miscellaneous routines.
00903     double LAPY2(const double x, const double y) const;
00904   };
00905 
00906   // END DOUBLE PARTIAL SPECIALIZATION DECLARATION //
00907 
00908   // BEGIN DOUBLE PARTIAL SPECIALIZATION IMPLEMENTATION //
00909 
00910 
00911   template<typename OrdinalType>
00912   void LAPACK<OrdinalType, double>::POTRF(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const
00913   {
00914     DPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00915   }
00916   
00917   template<typename OrdinalType>
00918   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
00919   {
00920     DPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00921   }
00922   
00923   template<typename OrdinalType>
00924   void LAPACK<OrdinalType, double>::POTRI(const char UPLO, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* info) const
00925   {
00926     DPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
00927   }
00928   
00929   template<typename OrdinalType>
00930     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
00931   {
00932     DPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00933   }
00934   
00935   template<typename OrdinalType>
00936   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
00937   {
00938     DPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
00939   }
00940   
00941   template<typename OrdinalType>
00942   void LAPACK<OrdinalType, double>::POEQU(const OrdinalType n, const double* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const
00943   {
00944     DPOEQU_F77(&n, A, &lda, S, scond, amax, info);
00945   }
00946   
00947   template<typename OrdinalType>
00948   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
00949   {
00950     DPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
00951   }
00952   
00953   template<typename OrdinalType>
00954     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 
00955   {
00956     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);
00957   }
00958   
00959   template<typename OrdinalType>
00960   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
00961   {
00962     DGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, &info);
00963   }
00964   
00965   template<typename OrdinalType>  
00966   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
00967   {
00968     DGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
00969   }
00970 
00971   template<typename OrdinalType>
00972   void LAPACK<OrdinalType,double>::GETRF(const OrdinalType m, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
00973   {
00974     DGETRF_F77(&m, &n, A, &lda, IPIV, info);
00975   }
00976   
00977   template<typename OrdinalType>
00978   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
00979   {
00980     DGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00981   }
00982   
00983   template<typename OrdinalType>
00984   void LAPACK<OrdinalType,double>::GETRI(const OrdinalType n, double* A, const OrdinalType lda, const OrdinalType* IPIV, double* WORK, const OrdinalType lwork, OrdinalType* info) const
00985   {
00986     DGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info);
00987   }
00988   
00989   template<typename OrdinalType>
00990   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
00991   {
00992     DGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, IWORK, info);
00993   }
00994   
00995   template<typename OrdinalType>
00996   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
00997   {
00998     DGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
00999   }
01000   
01001   template<typename OrdinalType>
01002   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
01003   {
01004     DGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
01005   }
01006   
01007   template<typename OrdinalType>
01008   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
01009   {
01010     DGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info);
01011   }
01012   
01013   template<typename OrdinalType>
01014   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
01015   {
01016     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);
01017   }
01018   
01019   template<typename OrdinalType>
01020   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
01021   {
01022     DSYTRD_F77(CHAR_MACRO(UPLO), &n, A, &lda, D, E, TAU, WORK, &lwork, info);
01023   }
01024 
01025   template<typename OrdinalType>
01026   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
01027   {
01028     DGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01029   }
01030 
01031   template<typename OrdinalType>
01032   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
01033   {
01034     DSPEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, AP, W, Z, &ldz, WORK, info);
01035   }
01036 
01037   template<typename OrdinalType>
01038   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
01039   {
01040     DSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info);
01041   }
01042 
01043   template<typename OrdinalType>
01044   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
01045   {
01046     DSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info);
01047   }
01048 
01049   template<typename OrdinalType>
01050   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
01051   {
01052     DSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01053   }
01054 
01055   template<typename OrdinalType>
01056   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
01057   {
01058     DHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, WR, WI, Z, &ldz, WORK, &lwork, info);
01059   }
01060   
01061   template<typename OrdinalType>
01062   void LAPACK<OrdinalType, double>::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, double* A, const OrdinalType lda, OrdinalType* sdim, double* WR, double* WI, double* VS, const OrdinalType ldvs, double* WORK, const OrdinalType lwork, OrdinalType* BWORK, OrdinalType* info) const    
01063   {
01064     DGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
01065   }
01066 
01067   template<typename OrdinalType>
01068   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
01069   {
01070     DGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, WORK, &lwork, info);
01071   }
01072 
01073   template<typename OrdinalType>
01074   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
01075   {
01076     DORMQR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
01077   }
01078 
01079   template<typename OrdinalType>
01080   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
01081   {
01082     DORGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info);
01083   }
01084   
01085   template<typename OrdinalType>
01086   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
01087   {
01088     DORGHR_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01089   }
01090   
01091   template<typename OrdinalType>
01092   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
01093   {
01094     DORMHR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &ilo, &ihi, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
01095   }
01096   
01097   template<typename OrdinalType>
01098   void LAPACK<OrdinalType, double>::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const double* T, const OrdinalType ldt, double* VL, const OrdinalType ldvl, double* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, double* WORK, OrdinalType* info) const
01099   {
01100     if(HOWMNY=='S')
01101       {
01102   *info = -3; // We do not support 'S' since it requires a logical function (yuck!)
01103       }
01104     else
01105       {
01106   DTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
01107       }
01108   }
01109   
01110   template<typename OrdinalType>
01111   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
01112   {
01113     DTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, WORK, info);
01114   }
01115 
01116   template<typename OrdinalType>
01117   double LAPACK<OrdinalType, double>::LARND( const OrdinalType idist, OrdinalType* seed ) const
01118   {
01119     return(DLARND_F77(&idist, seed));
01120   }
01121 
01122   template<typename OrdinalType>
01123   void LAPACK<OrdinalType, double>::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, double* v ) const
01124   {
01125     DLARNV_F77(&idist, seed, &n, v);
01126   }
01127 
01128   template<typename OrdinalType>
01129   double LAPACK<OrdinalType, double>::LAMCH(const char CMACH) const
01130   {
01131     return(DLAMCH_F77(CHAR_MACRO(CMACH)));
01132   }
01133 
01134   template<typename OrdinalType>
01135   int LAPACK<OrdinalType, double>::ILAENV( const int ispec, const string NAME, const string OPTS, const int N1, const int N2, const int N3, const int N4 ) const
01136   {
01137     unsigned int opts_length = OPTS.length();
01138     string temp_NAME = "d" + NAME;
01139     unsigned int name_length = temp_NAME.length();
01140 #if defined (INTEL_CXML)
01141     return ILAENV_F77(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4 );
01142 #else
01143     return ILAENV_F77(&ispec, &temp_NAME[0], &OPTS[0], &N1, &N2, &N3, &N4, name_length, opts_length );
01144 #endif
01145   }
01146  
01147   template<typename OrdinalType>
01148   double LAPACK<OrdinalType, double>::LAPY2(const double x, const double y) const
01149   {
01150     return DLAPY2_F77(&x, &y);
01151   }
01152 
01153   // END DOUBLE PARTIAL SPECIALIZATION IMPLEMENTATION //
01154 
01155 #ifdef HAVE_TEUCHOS_COMPLEX
01156 
01157   // BEGIN COMPLEX<FLOAT> PARTIAL SPECIALIZATION DECLARATION //
01158 
01159   template<typename OrdinalType>
01160   class LAPACK<OrdinalType, complex<float> >
01161   {    
01162   public:
01163     inline LAPACK(void) {};
01164     inline LAPACK(const LAPACK<OrdinalType, complex<float> >& LAPACK) {};
01165     inline virtual ~LAPACK(void) {};
01166 
01167     // Symmetric positive definite linear system routines
01168     void POTRF(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const;
01169     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;
01170     void POTRI(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const;
01171     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;
01172     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;
01173     void POEQU(const OrdinalType n, const complex<float>* A, const OrdinalType lda, float* S, float* scond, float* amax, OrdinalType* info) const;
01174     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;
01175     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; 
01176 
01177     // General Linear System Routines
01178     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;
01179     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;
01180     void GETRF(const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
01181     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;
01182     void GETRI(const OrdinalType n, complex<float>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<float>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01183     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;
01184     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;
01185     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;
01186     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;
01187     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;
01188     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;
01189 
01190     // Symmetric eigenvalue routines.
01191     void STEQR(const char COMPZ, const OrdinalType n, float* D, float* E, complex<float>* Z, const OrdinalType ldz, float* WORK, OrdinalType* info) const;
01192 
01193     // Hessenberg eigenvalue routines.
01194     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;
01195     void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* sdim, complex<float>* W, complex<float>* VS, const OrdinalType ldvs, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* BWORK, OrdinalType* info) const;    
01196     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;
01197 
01198     // Triangular matrix routines.
01199     void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<float>* T, const OrdinalType ldt, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<float>* WORK, float* RWORK, OrdinalType* info) const;
01200     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;
01201 
01202     // Random number generators
01203     complex<float> LARND( const OrdinalType idist, OrdinalType* seed ) const;
01204     void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<float>* v ) const;    
01205 
01206     // Machine characteristics
01207     int ILAENV( const int ispec, const string NAME, const string OPTS, const int N1 = -1, const int N2 = -1, const int N3 = -1, const int N4 = -1 ) const;
01208 
01209   };
01210 
01211   // END COMPLEX<FLOAT> PARTIAL SPECIALIZATION DECLARATION //
01212 
01213   // BEGIN COMPLEX<FLOAT> PARTIAL SPECIALIZATION IMPLEMENTATION //
01214 
01215   template<typename OrdinalType>
01216   void LAPACK<OrdinalType, complex<float> >::POTRF(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const
01217   {
01218     CPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01219   }
01220   
01221   template<typename OrdinalType>
01222   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
01223   {
01224     CPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01225   }
01226   
01227   template<typename OrdinalType>
01228   void LAPACK<OrdinalType, complex<float> >::POTRI(const char UPLO, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* info) const
01229   {
01230     CPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01231   }
01232   
01233   template<typename OrdinalType>
01234   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
01235   {
01236     CPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01237   }
01238   
01239   template<typename OrdinalType>
01240   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
01241   {
01242     CPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01243   }
01244   
01245   template<typename OrdinalType>
01246   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
01247   {
01248     CPOEQU_F77(&n, A, &lda, S, scond, amax, info);
01249   }
01250   
01251   template<typename OrdinalType>
01252   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
01253   {
01254     CPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01255   }
01256   
01257   template<typename OrdinalType>
01258   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
01259   {
01260     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);
01261   }
01262   
01263   template<typename OrdinalType>
01264   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
01265   {
01266     CGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
01267   }
01268   
01269   template<typename OrdinalType>  
01270   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
01271   {
01272     CGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
01273   }
01274 
01275   template<typename OrdinalType>
01276   void LAPACK<OrdinalType,complex<float> >::GETRF(const OrdinalType m, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
01277   {
01278     CGETRF_F77(&m, &n, A, &lda, IPIV, info);
01279   }
01280   
01281   template<typename OrdinalType>
01282   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
01283   {
01284     CGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01285   }
01286   
01287   template<typename OrdinalType>
01288   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
01289   {
01290     CGETRI_F77(&n, A, &lda, IPIV, WORK, lwork, info);
01291   }
01292   
01293   template<typename OrdinalType>
01294   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
01295   {
01296     CGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01297   }
01298   
01299   template<typename OrdinalType>
01300   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
01301   {
01302     CGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01303   }
01304   
01305   template<typename OrdinalType>
01306   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
01307   {
01308     CGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
01309   }
01310   
01311   template<typename OrdinalType>
01312   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
01313   {
01314     CGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01315   }
01316   
01317   template<typename OrdinalType>
01318   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
01319   {
01320     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);
01321   }
01322   
01323   template<typename OrdinalType>
01324   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
01325   {
01326     CGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01327   }
01328   
01329   template<typename OrdinalType>
01330   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
01331   {
01332     CSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01333   }
01334 
01335   template<typename OrdinalType>
01336   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
01337   {
01338     CHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
01339   }
01340 
01341   template<typename OrdinalType>
01342   void LAPACK<OrdinalType, complex<float> >::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<float>* A, const OrdinalType lda, OrdinalType* sdim, complex<float>* W, complex<float>* VS, const OrdinalType ldvs, complex<float>* WORK, const OrdinalType lwork, float* RWORK, OrdinalType* BWORK, OrdinalType* info) const    
01343   {
01344     CGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
01345   }
01346 
01347   template<typename OrdinalType>
01348   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
01349   {
01350     CGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
01351   }
01352 
01353   template<typename OrdinalType>
01354   void LAPACK<OrdinalType, complex<float> >::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<float>* T, const OrdinalType ldt, complex<float>* VL, const OrdinalType ldvl, complex<float>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<float>* WORK, float* RWORK, OrdinalType* info) const
01355   {
01356     if(HOWMNY=='S')
01357       {
01358   *info = -3; // We do not support 'S' since it requires a logical function (yuck!)
01359       }
01360     else
01361       {
01362   CTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
01363       }
01364   }
01365   
01366   template<typename OrdinalType>
01367   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
01368   {
01369     CTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, info);
01370   }
01371 
01372   template<typename OrdinalType>
01373   complex<float> LAPACK<OrdinalType, complex<float> >::LARND( const OrdinalType idist, OrdinalType* seed ) const
01374   {
01375     return(CLARND_F77(&idist, seed));
01376   }
01377 
01378   template<typename OrdinalType>
01379   void LAPACK<OrdinalType, complex<float> >::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<float>* v ) const
01380   {
01381     CLARNV_F77(&idist, seed, &n, v);
01382   }
01383 
01384   template<typename OrdinalType>
01385   int LAPACK<OrdinalType, complex<float> >::ILAENV( const int ispec, const string NAME, const string OPTS, const int N1, const int N2, const int N3, const int N4 ) const
01386   {
01387     unsigned int opts_length = OPTS.length();
01388     string temp_NAME = "c" + NAME;
01389     unsigned int name_length = temp_NAME.length();
01390 #if defined (INTEL_CXML)
01391     return ILAENV_F77(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4 );
01392 #else
01393     return ILAENV_F77(&ispec, &temp_NAME[0], &OPTS[0], &N1, &N2, &N3, &N4, name_length, opts_length );
01394 #endif
01395   }
01396 
01397   // END COMPLEX<FLOAT> PARTIAL SPECIALIZATION IMPLEMENTATION //
01398 
01399   // BEGIN COMPLEX<DOUBLE> PARTIAL SPECIALIZATION DECLARATION //
01400 
01401   template<typename OrdinalType>
01402   class LAPACK<OrdinalType, complex<double> >
01403   {    
01404   public:
01405     inline LAPACK(void) {};
01406     inline LAPACK(const LAPACK<OrdinalType, complex<double> >& LAPACK) {};
01407     inline virtual ~LAPACK(void) {};
01408 
01409     // Symmetric positive definite linear system routines
01410     void POTRF(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const;
01411     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;
01412     void POTRI(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const;
01413     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;
01414     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;
01415     void POEQU(const OrdinalType n, const complex<double>* A, const OrdinalType lda, double* S, double* scond, double* amax, OrdinalType* info) const;
01416     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;
01417     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; 
01418 
01419     // General Linear System Routines
01420     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;
01421     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;
01422     void GETRF(const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const;
01423     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;
01424     void GETRI(const OrdinalType n, complex<double>* A, const OrdinalType lda, const OrdinalType* IPIV, complex<double>* WORK, const OrdinalType lwork, OrdinalType* info) const;
01425     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;
01426     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;
01427     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;
01428     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;
01429     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;
01430     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;
01431 
01432     // Symmetric eigenvalue routines.
01433     void STEQR(const char COMPZ, const OrdinalType n, double* D, double* E, complex<double>* Z, const OrdinalType ldz, double* WORK, OrdinalType* info) const;
01434 
01435     // Hessenberg eigenvalue routines.
01436     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;
01437     void GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* sdim, complex<double>* W, complex<double>* VS, const OrdinalType ldvs, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* BWORK, OrdinalType* info) const;    
01438     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;
01439 
01440     // Triangular matrix routines.
01441     void TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<double>* T, const OrdinalType ldt, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<double>* WORK, double* RWORK, OrdinalType* info) const;
01442     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;
01443 
01444     // Random number generators
01445     complex<double> LARND( const OrdinalType idist, OrdinalType* seed ) const;
01446     void LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<double>* v ) const;    
01447 
01448     // Machine characteristics
01449     int ILAENV( const int ispec, const string NAME, const string OPTS, const int N1 = -1, const int N2 = -1, const int N3 = -1, const int N4 = -1 ) const;
01450 
01451   };
01452 
01453   // END COMPLEX<DOUBLE> PARTIAL SPECIALIZATION DECLARATION //
01454 
01455   // BEGIN COMPLEX<DOUBLE> PARTIAL SPECIALIZATION IMPLEMENTATION //
01456 
01457   template<typename OrdinalType>
01458   void LAPACK<OrdinalType, complex<double> >::POTRF(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const
01459   {
01460     ZPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01461   }
01462   
01463   template<typename OrdinalType>
01464   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
01465   {
01466     ZPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01467   }
01468   
01469   template<typename OrdinalType>
01470   void LAPACK<OrdinalType, complex<double> >::POTRI(const char UPLO, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* info) const
01471   {
01472     ZPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
01473   }
01474   
01475   template<typename OrdinalType>
01476   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
01477   {
01478     ZPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01479   }
01480   
01481   template<typename OrdinalType>
01482   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
01483   {
01484     ZPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
01485   }
01486   
01487   template<typename OrdinalType>
01488   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
01489   {
01490     ZPOEQU_F77(&n, A, &lda, S, scond, amax, info);
01491   }
01492   
01493   template<typename OrdinalType>
01494   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
01495   {
01496     ZPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01497   }
01498   
01499   template<typename OrdinalType>
01500   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
01501   {
01502     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);
01503   }
01504   
01505   template<typename OrdinalType>
01506   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
01507   {
01508     ZGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
01509   }
01510   
01511   template<typename OrdinalType>  
01512   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
01513   {
01514     ZGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
01515   }
01516 
01517   template<typename OrdinalType>
01518   void LAPACK<OrdinalType,complex<double> >::GETRF(const OrdinalType m, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* IPIV, OrdinalType* info) const
01519   {
01520     ZGETRF_F77(&m, &n, A, &lda, IPIV, info);
01521   }
01522   
01523   template<typename OrdinalType>
01524   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
01525   {
01526     ZGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01527   }
01528   
01529   template<typename OrdinalType>
01530   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
01531   {
01532     ZGETRI_F77(&n, A, &lda, IPIV, WORK, lwork, info);
01533   }
01534   
01535   template<typename OrdinalType>
01536   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
01537   {
01538     ZGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
01539   }
01540   
01541   template<typename OrdinalType>
01542   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
01543   {
01544     ZGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
01545   }
01546   
01547   template<typename OrdinalType>
01548   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
01549   {
01550     ZGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
01551   }
01552   
01553   template<typename OrdinalType>
01554   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
01555   {
01556     ZGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
01557   }
01558   
01559   template<typename OrdinalType>
01560   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
01561   {
01562     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);
01563   }
01564   
01565   template<typename OrdinalType>
01566   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
01567   {
01568     ZGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
01569   }
01570   
01571   template<typename OrdinalType>
01572   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
01573   {
01574     ZSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
01575   }
01576 
01577   template<typename OrdinalType>
01578   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
01579   {
01580     ZHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
01581   }
01582 
01583   template<typename OrdinalType>
01584   void LAPACK<OrdinalType, complex<double> >::GEES(const char JOBVS, const char SORT, OrdinalType* SELECT, const OrdinalType n, complex<double>* A, const OrdinalType lda, OrdinalType* sdim, complex<double>* W, complex<double>* VS, const OrdinalType ldvs, complex<double>* WORK, const OrdinalType lwork, double* RWORK, OrdinalType* BWORK, OrdinalType* info) const    
01585   {
01586     ZGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), SELECT, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
01587   }
01588 
01589   template<typename OrdinalType>
01590   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
01591   {
01592     ZGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
01593   }
01594 
01595   template<typename OrdinalType>
01596   void LAPACK<OrdinalType, complex<double> >::TREVC(const char SIDE, const char HOWMNY, OrdinalType* SELECT, const OrdinalType n, const complex<double>* T, const OrdinalType ldt, complex<double>* VL, const OrdinalType ldvl, complex<double>* VR, const OrdinalType ldvr, const OrdinalType mm, OrdinalType* m, complex<double>* WORK, double* RWORK, OrdinalType* info) const
01597   {
01598     if(HOWMNY=='S')
01599       {
01600   *info = -3; // We do not support 'S' since it requires a logical function (yuck!)
01601       }
01602     else
01603       {
01604   ZTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
01605       }
01606   }
01607   
01608   template<typename OrdinalType>
01609   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
01610   {
01611     ZTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, &ifst, &ilst, info);
01612   }
01613 
01614   template<typename OrdinalType>
01615   complex<double> LAPACK<OrdinalType, complex<double> >::LARND( const OrdinalType idist, OrdinalType* seed ) const
01616   {
01617     return(ZLARND_F77(&idist, seed));
01618   }
01619 
01620   template<typename OrdinalType>
01621   void LAPACK<OrdinalType, complex<double> >::LARNV( const OrdinalType idist, OrdinalType* seed, const OrdinalType n, complex<double>* v ) const
01622   {
01623     ZLARNV_F77(&idist, seed, &n, v);
01624   }
01625 
01626   template<typename OrdinalType>
01627   int LAPACK<OrdinalType, complex<double> >::ILAENV( const int ispec, const string NAME, const string OPTS, const int N1, const int N2, const int N3, const int N4 ) const
01628   {
01629     unsigned int opts_length = OPTS.length();
01630     string temp_NAME = "z" + NAME;
01631     unsigned int name_length = temp_NAME.length();
01632 #if defined (INTEL_CXML)
01633     return ILAENV_F77(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4 );
01634 #else
01635     return ILAENV_F77(&ispec, &temp_NAME[0], &OPTS[0], &N1, &N2, &N3, &N4, name_length, opts_length );
01636 #endif
01637   }
01638 
01639   // END COMPLEX<DOUBLE> PARTIAL SPECIALIZATION IMPLEMENTATION //
01640 
01641 #endif // HAVE_TEUCHOS_COMPLEX
01642 
01643 #endif // DOXYGEN_SHOULD_SKIP_THIS
01644 
01645 } // namespace Teuchos
01646 
01647 #endif // _TEUCHOS_LAPACK_HPP_

Generated on Thu Sep 18 12:41:17 2008 for Teuchos - Trilinos Tools Package by doxygen 1.3.9.1