[2322] | 1 | #include <iostream>
|
---|
[2567] | 2 | #include <math.h>
|
---|
[775] | 3 | #include "intflapack.h"
|
---|
[1342] | 4 | #include "tvector.h"
|
---|
| 5 | #include "tmatrix.h"
|
---|
[814] | 6 | #include <typeinfo>
|
---|
[775] | 7 |
|
---|
[2567] | 8 | #define GARDMEM 5
|
---|
| 9 |
|
---|
[2556] | 10 | /*************** Pour memoire (Christophe) ***************
|
---|
| 11 | Les dispositions memoires (FORTRAN) pour les vecteurs et matrices LAPACK:
|
---|
| 12 |
|
---|
| 13 | 1./ --- REAL X(N):
|
---|
| 14 | if an array X of dimension (N) holds a vector x,
|
---|
| 15 | then X(i) holds "x_i" for i=1,...,N
|
---|
| 16 |
|
---|
| 17 | 2./ --- REAL A(LDA,N):
|
---|
| 18 | if a two-dimensional array A of dimension (LDA,N) holds an m-by-n matrix A,
|
---|
| 19 | then A(i,j) holds "a_ij" for i=1,...,m et j=1,...,n (LDA must be at least m).
|
---|
| 20 | Note that array arguments are usually declared in the software as assumed-size
|
---|
| 21 | arrays (last dimension *), for example: REAL A(LDA,*)
|
---|
| 22 | --- Rangement en memoire:
|
---|
| 23 | | 11 12 13 14 |
|
---|
| 24 | Ex: Real A(4,4): A = | 21 22 23 24 |
|
---|
| 25 | | 31 32 33 34 |
|
---|
| 26 | | 41 42 43 44 |
|
---|
| 27 | memoire: {11 21 31 41} {12 22 32 42} {13 23 33 43} {14 24 34 44}
|
---|
| 28 | First indice (line) "i" varies then the second (column):
|
---|
| 29 | (put all the first column, then put all the second column,
|
---|
| 30 | ..., then put all the last column)
|
---|
| 31 | ***********************************************************/
|
---|
| 32 |
|
---|
[1424] | 33 | /*!
|
---|
| 34 | \defgroup LinAlg LinAlg module
|
---|
| 35 | This module contains classes and functions for complex linear
|
---|
| 36 | algebra on arrays. This module is intended mainly to have
|
---|
| 37 | classes implementing C++ interfaces between Sophya objects
|
---|
| 38 | and external linear algebra libraries, such as LAPACK.
|
---|
| 39 | */
|
---|
| 40 |
|
---|
| 41 | /*!
|
---|
| 42 | \class SOPHYA::LapackServer
|
---|
| 43 | \ingroup LinAlg
|
---|
| 44 | This class implements an interface to LAPACK library driver routines.
|
---|
| 45 | The LAPACK (Linear Algebra PACKage) is a collection high performance
|
---|
| 46 | routines to solve common problems in numerical linear algebra.
|
---|
| 47 | its is available from http://www.netlib.org.
|
---|
| 48 |
|
---|
| 49 | The present version of our LapackServer (Feb 2001) provides only
|
---|
| 50 | interfaces for the linear system solver and singular value
|
---|
| 51 | decomposition (SVD). Only arrays with BaseArray::FortranMemoryMapping
|
---|
| 52 | can be handled by LapackServer. LapackServer can be instanciated
|
---|
| 53 | for simple and double precision real or complex array types.
|
---|
| 54 |
|
---|
| 55 | The example below shows solving a linear system A*X = B
|
---|
| 56 |
|
---|
| 57 | \code
|
---|
| 58 | #include "intflapack.h"
|
---|
| 59 | // ...
|
---|
| 60 | // Use FortranMemoryMapping as default
|
---|
| 61 | BaseArray::SetDefaultMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
| 62 | // Create an fill the arrays A and B
|
---|
| 63 | int n = 20;
|
---|
| 64 | Matrix A(n, n);
|
---|
| 65 | A = RandomSequence();
|
---|
| 66 | Vector X(n),B(n);
|
---|
| 67 | X = RandomSequence();
|
---|
| 68 | B = A*X;
|
---|
| 69 | // Solve the linear system A*X = B
|
---|
| 70 | LapackServer<r_8> lps;
|
---|
| 71 | lps.LinSolve(A,B);
|
---|
| 72 | // We get the result in B, which should be equal to X ...
|
---|
| 73 | // Compute the difference B-X ;
|
---|
| 74 | Vector diff = B-X;
|
---|
| 75 | \endcode
|
---|
| 76 |
|
---|
| 77 | */
|
---|
| 78 |
|
---|
[2556] | 79 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[775] | 80 | extern "C" {
|
---|
[2554] | 81 | // Le calculateur de workingspace
|
---|
| 82 | int_4 ilaenv_(int_4 *ispec,char *name,char *opts,int_4 *n1,int_4 *n2,int_4 *n3,int_4 *n4,
|
---|
| 83 | int_4 nc1,int_4 nc2);
|
---|
| 84 |
|
---|
[1342] | 85 | // Drivers pour resolution de systemes lineaires
|
---|
| 86 | void sgesv_(int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
| 87 | int_4* ipiv, r_4* b, int_4* ldb, int_4* info);
|
---|
| 88 | void dgesv_(int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
| 89 | int_4* ipiv, r_8* b, int_4* ldb, int_4* info);
|
---|
| 90 | void cgesv_(int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
| 91 | int_4* ipiv, complex<r_4>* b, int_4* ldb, int_4* info);
|
---|
| 92 | void zgesv_(int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
| 93 | int_4* ipiv, complex<r_8>* b, int_4* ldb, int_4* info);
|
---|
| 94 |
|
---|
[2554] | 95 | // Drivers pour resolution de systemes lineaires symetriques
|
---|
| 96 | void ssysv_(char* uplo, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
| 97 | int_4* ipiv, r_4* b, int_4* ldb,
|
---|
| 98 | r_4* work, int_4* lwork, int_4* info);
|
---|
| 99 | void dsysv_(char* uplo, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
| 100 | int_4* ipiv, r_8* b, int_4* ldb,
|
---|
| 101 | r_8* work, int_4* lwork, int_4* info);
|
---|
| 102 | void csysv_(char* uplo, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
| 103 | int_4* ipiv, complex<r_4>* b, int_4* ldb,
|
---|
| 104 | complex<r_4>* work, int_4* lwork, int_4* info);
|
---|
| 105 | void zsysv_(char* uplo, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
| 106 | int_4* ipiv, complex<r_8>* b, int_4* ldb,
|
---|
| 107 | complex<r_8>* work, int_4* lwork, int_4* info);
|
---|
| 108 |
|
---|
| 109 | // Driver pour resolution de systemes au sens de Xi2
|
---|
[1494] | 110 | void sgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
| 111 | r_4* b, int_4* ldb, r_4* work, int_4* lwork, int_4* info);
|
---|
| 112 | void dgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
| 113 | r_8* b, int_4* ldb, r_8* work, int_4* lwork, int_4* info);
|
---|
| 114 | void cgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
| 115 | complex<r_4>* b, int_4* ldb, complex<r_4>* work, int_4* lwork, int_4* info);
|
---|
| 116 | void zgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
| 117 | complex<r_8>* b, int_4* ldb, complex<r_8>* work, int_4* lwork, int_4* info);
|
---|
| 118 |
|
---|
[2567] | 119 | // Driver pour resolution de systemes au sens de Xi2 par SVD Divide & Conquer
|
---|
| 120 | void sgelsd_(int_4* m,int_4* n,int_4* nrhs,r_4* a,int_4* lda,
|
---|
| 121 | r_4* b,int_4* ldb,r_4* s,r_4* rcond,int_4* rank,
|
---|
| 122 | r_4* work,int_4* lwork,int_4* iwork,int_4* info);
|
---|
| 123 | void dgelsd_(int_4* m,int_4* n,int_4* nrhs,r_8* a,int_4* lda,
|
---|
| 124 | r_8* b,int_4* ldb,r_8* s,r_8* rcond,int_4* rank,
|
---|
| 125 | r_8* work,int_4* lwork,int_4* iwork,int_4* info);
|
---|
| 126 | void cgelsd_(int_4* m,int_4* n,int_4* nrhs,complex<r_4>* a,int_4* lda,
|
---|
| 127 | complex<r_4>* b,int_4* ldb,r_4* s,r_4* rcond,int_4* rank,
|
---|
| 128 | complex<r_4>* work,int_4* lwork,r_4* rwork,int_4* iwork,int_4* info);
|
---|
| 129 | void zgelsd_(int_4* m,int_4* n,int_4* nrhs,complex<r_8>* a,int_4* lda,
|
---|
| 130 | complex<r_8>* b,int_4* ldb,r_8* s,r_8* rcond,int_4* rank,
|
---|
| 131 | complex<r_8>* work,int_4* lwork,r_8* rwork,int_4* iwork,int_4* info);
|
---|
| 132 |
|
---|
[1342] | 133 | // Driver pour decomposition SVD
|
---|
| 134 | void sgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_4* a, int_4* lda,
|
---|
| 135 | r_4* s, r_4* u, int_4* ldu, r_4* vt, int_4* ldvt,
|
---|
| 136 | r_4* work, int_4* lwork, int_4* info);
|
---|
| 137 | void dgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_8* a, int_4* lda,
|
---|
| 138 | r_8* s, r_8* u, int_4* ldu, r_8* vt, int_4* ldvt,
|
---|
| 139 | r_8* work, int_4* lwork, int_4* info);
|
---|
| 140 | void cgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_4>* a, int_4* lda,
|
---|
[2559] | 141 | r_4* s, complex<r_4>* u, int_4* ldu, complex<r_4>* vt, int_4* ldvt,
|
---|
| 142 | complex<r_4>* work, int_4* lwork, r_4* rwork, int_4* info);
|
---|
[1342] | 143 | void zgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_8>* a, int_4* lda,
|
---|
[2559] | 144 | r_8* s, complex<r_8>* u, int_4* ldu, complex<r_8>* vt, int_4* ldvt,
|
---|
| 145 | complex<r_8>* work, int_4* lwork, r_8* rwork, int_4* info);
|
---|
[2556] | 146 |
|
---|
[2561] | 147 | // Driver pour decomposition SVD Divide and Conquer
|
---|
| 148 | void sgesdd_(char* jobz, int_4* m, int_4* n, r_4* a, int_4* lda,
|
---|
| 149 | r_4* s, r_4* u, int_4* ldu, r_4* vt, int_4* ldvt,
|
---|
| 150 | r_4* work, int_4* lwork, int_4* iwork, int_4* info);
|
---|
| 151 | void dgesdd_(char* jobz, int_4* m, int_4* n, r_8* a, int_4* lda,
|
---|
| 152 | r_8* s, r_8* u, int_4* ldu, r_8* vt, int_4* ldvt,
|
---|
| 153 | r_8* work, int_4* lwork, int_4* iwork, int_4* info);
|
---|
| 154 | void cgesdd_(char* jobz, int_4* m, int_4* n, complex<r_4>* a, int_4* lda,
|
---|
| 155 | r_4* s, complex<r_4>* u, int_4* ldu, complex<r_4>* vt, int_4* ldvt,
|
---|
| 156 | complex<r_4>* work, int_4* lwork, r_4* rwork, int_4* iwork, int_4* info);
|
---|
| 157 | void zgesdd_(char* jobz, int_4* m, int_4* n, complex<r_8>* a, int_4* lda,
|
---|
| 158 | r_8* s, complex<r_8>* u, int_4* ldu, complex<r_8>* vt, int_4* ldvt,
|
---|
| 159 | complex<r_8>* work, int_4* lwork, r_8* rwork, int_4* iwork, int_4* info);
|
---|
| 160 |
|
---|
[2556] | 161 | // Driver pour eigen decomposition for symetric/hermitian matrices
|
---|
| 162 | void ssyev_(char* jobz, char* uplo, int_4* n, r_4* a, int_4* lda, r_4* w,
|
---|
| 163 | r_4* work, int_4 *lwork, int_4* info);
|
---|
| 164 | void dsyev_(char* jobz, char* uplo, int_4* n, r_8* a, int_4* lda, r_8* w,
|
---|
| 165 | r_8* work, int_4 *lwork, int_4* info);
|
---|
| 166 | void cheev_(char* jobz, char* uplo, int_4* n, complex<r_4>* a, int_4* lda, r_4* w,
|
---|
| 167 | complex<r_4>* work, int_4 *lwork, r_4* rwork, int_4* info);
|
---|
| 168 | void zheev_(char* jobz, char* uplo, int_4* n, complex<r_8>* a, int_4* lda, r_8* w,
|
---|
| 169 | complex<r_8>* work, int_4 *lwork, r_8* rwork, int_4* info);
|
---|
| 170 |
|
---|
| 171 | // Driver pour eigen decomposition for general squared matrices
|
---|
| 172 | void sgeev_(char* jobl, char* jobvr, int_4* n, r_4* a, int_4* lda, r_4* wr, r_4* wi,
|
---|
| 173 | r_4* vl, int_4* ldvl, r_4* vr, int_4* ldvr,
|
---|
| 174 | r_4* work, int_4 *lwork, int_4* info);
|
---|
| 175 | void dgeev_(char* jobl, char* jobvr, int_4* n, r_8* a, int_4* lda, r_8* wr, r_8* wi,
|
---|
| 176 | r_8* vl, int_4* ldvl, r_8* vr, int_4* ldvr,
|
---|
| 177 | r_8* work, int_4 *lwork, int_4* info);
|
---|
| 178 | void cgeev_(char* jobl, char* jobvr, int_4* n, complex<r_4>* a, int_4* lda, complex<r_4>* w,
|
---|
| 179 | complex<r_4>* vl, int_4* ldvl, complex<r_4>* vr, int_4* ldvr,
|
---|
| 180 | complex<r_4>* work, int_4 *lwork, r_4* rwork, int_4* info);
|
---|
| 181 | void zgeev_(char* jobl, char* jobvr, int_4* n, complex<r_8>* a, int_4* lda, complex<r_8>* w,
|
---|
| 182 | complex<r_8>* vl, int_4* ldvl, complex<r_8>* vr, int_4* ldvr,
|
---|
| 183 | complex<r_8>* work, int_4 *lwork, r_8* rwork, int_4* info);
|
---|
| 184 |
|
---|
[775] | 185 | }
|
---|
| 186 |
|
---|
[1342] | 187 | // -------------- Classe LapackServer<T> --------------
|
---|
| 188 |
|
---|
[2556] | 189 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[814] | 190 | template <class T>
|
---|
[1344] | 191 | LapackServer<T>::LapackServer()
|
---|
[1342] | 192 | {
|
---|
| 193 | SetWorkSpaceSizeFactor();
|
---|
| 194 | }
|
---|
| 195 |
|
---|
| 196 | template <class T>
|
---|
[1344] | 197 | LapackServer<T>::~LapackServer()
|
---|
[1342] | 198 | {
|
---|
| 199 | }
|
---|
| 200 |
|
---|
[2556] | 201 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[2554] | 202 | template <class T>
|
---|
| 203 | int_4 LapackServer<T>::ilaenv_en_C(int_4 ispec,char *name,char *opts,int_4 n1,int_4 n2,int_4 n3,int_4 n4)
|
---|
| 204 | {
|
---|
[2572] | 205 | int_4 nc1 = strlen(name), nc2 = strlen(opts), rc=0;
|
---|
[2554] | 206 | rc = ilaenv_(&ispec,name,opts,&n1,&n2,&n3,&n4,nc1,nc2);
|
---|
| 207 | //cout<<"ilaenv_en_C("<<ispec<<","<<name<<"("<<nc1<<"),"<<opts<<"("<<nc2<<"),"
|
---|
| 208 | // <<n1<<","<<n2<<","<<n3<<","<<n4<<") = "<<rc<<endl;
|
---|
| 209 | return rc;
|
---|
| 210 | }
|
---|
| 211 |
|
---|
[2572] | 212 | template <class T>
|
---|
| 213 | int_4 LapackServer<T>::type2i4(void *val,int nbytes)
|
---|
| 214 | // Retourne un entier contenant la valeur contenue dans val
|
---|
| 215 | // - nbytes = nombre de bytes dans le contenu de val
|
---|
| 216 | // ex: r_4 x = 3.4; type2i4(&x,4) -> 3
|
---|
| 217 | // ex: r_8 x = 3.4; type2i4(&x,8) -> 3
|
---|
| 218 | // ex: complex<r_4> x(3.4,7.8); type2i4(&x,4) -> 3
|
---|
| 219 | // ex: complex<r_8> x(3.4,7.8); type2i4(&x,8) -> 3
|
---|
| 220 | {
|
---|
| 221 | r_4* x4; r_8* x8; int_4 lw=0;
|
---|
| 222 | if(nbytes==4) {x4 = (r_4*)val; lw = (int_4)(*x4);}
|
---|
| 223 | else {x8 = (r_8*)val; lw = (int_4)(*x8);}
|
---|
| 224 | return lw;
|
---|
| 225 | }
|
---|
| 226 |
|
---|
[2556] | 227 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[2563] | 228 | //! Interface to Lapack linear system solver driver s/d/c/zgesv().
|
---|
| 229 | /*! Solve the linear system a * x = b using LU factorization.
|
---|
| 230 | Input arrays should have FortranMemory mapping (column packed).
|
---|
[1424] | 231 | \param a : input matrix, overwritten on output
|
---|
| 232 | \param b : input-output, input vector b, contains x on exit
|
---|
| 233 | \return : return code from lapack driver _gesv()
|
---|
| 234 | */
|
---|
[1342] | 235 | template <class T>
|
---|
[1042] | 236 | int LapackServer<T>::LinSolve(TArray<T>& a, TArray<T> & b)
|
---|
[814] | 237 | {
|
---|
| 238 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
| 239 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b NbDimensions() != 2"));
|
---|
| 240 |
|
---|
[1342] | 241 | int_4 rowa = a.RowsKA();
|
---|
| 242 | int_4 cola = a.ColsKA();
|
---|
| 243 | int_4 rowb = b.RowsKA();
|
---|
| 244 | int_4 colb = b.ColsKA();
|
---|
[814] | 245 | if ( a.Size(rowa) != a.Size(cola))
|
---|
| 246 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Not a square Array"));
|
---|
[1042] | 247 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
[814] | 248 | throw(SzMismatchError("LapackServer::LinSolve(a,b) RowSize(a <> b) "));
|
---|
| 249 |
|
---|
| 250 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
[1342] | 251 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b Not Column Packed"));
|
---|
[814] | 252 |
|
---|
| 253 | int_4 n = a.Size(rowa);
|
---|
| 254 | int_4 nrhs = b.Size(colb);
|
---|
| 255 | int_4 lda = a.Step(cola);
|
---|
| 256 | int_4 ldb = b.Step(colb);
|
---|
| 257 | int_4 info;
|
---|
| 258 | int_4* ipiv = new int_4[n];
|
---|
| 259 |
|
---|
| 260 | if (typeid(T) == typeid(r_4) )
|
---|
| 261 | sgesv_(&n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb, &info);
|
---|
| 262 | else if (typeid(T) == typeid(r_8) )
|
---|
| 263 | dgesv_(&n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb, &info);
|
---|
| 264 | else if (typeid(T) == typeid(complex<r_4>) )
|
---|
| 265 | cgesv_(&n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
---|
| 266 | (complex<r_4> *)b.Data(), &ldb, &info);
|
---|
| 267 | else if (typeid(T) == typeid(complex<r_8>) )
|
---|
| 268 | zgesv_(&n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
---|
| 269 | (complex<r_8> *)b.Data(), &ldb, &info);
|
---|
| 270 | else {
|
---|
| 271 | delete[] ipiv;
|
---|
| 272 | string tn = typeid(T).name();
|
---|
| 273 | cerr << " LapackServer::LinSolve(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 274 | throw TypeMismatchExc("LapackServer::LinSolve(a,b) - Unsupported DataType (T)");
|
---|
| 275 | }
|
---|
| 276 | delete[] ipiv;
|
---|
[1042] | 277 | return(info);
|
---|
[814] | 278 | }
|
---|
| 279 |
|
---|
[2556] | 280 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[2563] | 281 | //! Interface to Lapack linear system solver driver s/d/c/zsysv().
|
---|
| 282 | /*! Solve the linear system a * x = b with a symetric matrix using LU factorization.
|
---|
| 283 | Input arrays should have FortranMemory mapping (column packed).
|
---|
[2554] | 284 | \param a : input matrix symetric , overwritten on output
|
---|
| 285 | \param b : input-output, input vector b, contains x on exit
|
---|
[2563] | 286 | \return : return code from lapack driver _sysv()
|
---|
[2554] | 287 | */
|
---|
| 288 | template <class T>
|
---|
| 289 | int LapackServer<T>::LinSolveSym(TArray<T>& a, TArray<T> & b)
|
---|
| 290 | // --- REMARQUES DE CMV ---
|
---|
| 291 | // 1./ contrairement a ce qui est dit dans la doc, il s'agit
|
---|
| 292 | // de matrices SYMETRIQUES complexes et non HERMITIENNES !!!
|
---|
| 293 | // 2./ pourquoi les routines de LinSolve pour des matrices symetriques
|
---|
[2556] | 294 | // sont plus de deux fois plus lentes que les LinSolve generales sur OSF
|
---|
| 295 | // et sensiblement plus lentes sous Linux ???
|
---|
[2554] | 296 | {
|
---|
| 297 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
| 298 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Or b NbDimensions() != 2"));
|
---|
| 299 | int_4 rowa = a.RowsKA();
|
---|
| 300 | int_4 cola = a.ColsKA();
|
---|
| 301 | int_4 rowb = b.RowsKA();
|
---|
| 302 | int_4 colb = b.ColsKA();
|
---|
| 303 | if ( a.Size(rowa) != a.Size(cola))
|
---|
| 304 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Not a square Array"));
|
---|
| 305 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
| 306 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) RowSize(a <> b) "));
|
---|
| 307 |
|
---|
| 308 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
| 309 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Or b Not Column Packed"));
|
---|
| 310 |
|
---|
| 311 | int_4 n = a.Size(rowa);
|
---|
| 312 | int_4 nrhs = b.Size(colb);
|
---|
| 313 | int_4 lda = a.Step(cola);
|
---|
| 314 | int_4 ldb = b.Step(colb);
|
---|
| 315 | int_4 info = 0;
|
---|
| 316 | int_4* ipiv = new int_4[n];
|
---|
| 317 | int_4 lwork = -1;
|
---|
| 318 | T * work = NULL;
|
---|
[2572] | 319 | T wkget[2];
|
---|
[2554] | 320 |
|
---|
| 321 | char uplo = 'U'; // char uplo = 'L';
|
---|
| 322 | char struplo[5]; struplo[0] = uplo; struplo[1] = '\0';
|
---|
| 323 |
|
---|
| 324 | if (typeid(T) == typeid(r_4) ) {
|
---|
| 325 | ssysv_(&uplo, &n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb,
|
---|
[2572] | 326 | (r_4 *)wkget, &lwork, &info);
|
---|
| 327 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 328 | ssysv_(&uplo, &n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb,
|
---|
[2554] | 329 | (r_4 *)work, &lwork, &info);
|
---|
| 330 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
| 331 | dsysv_(&uplo, &n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb,
|
---|
[2572] | 332 | (r_8 *)wkget, &lwork, &info);
|
---|
| 333 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 334 | dsysv_(&uplo, &n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb,
|
---|
[2554] | 335 | (r_8 *)work, &lwork, &info);
|
---|
| 336 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
| 337 | csysv_(&uplo, &n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
---|
| 338 | (complex<r_4> *)b.Data(), &ldb,
|
---|
[2572] | 339 | (complex<r_4> *)wkget, &lwork, &info);
|
---|
| 340 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 341 | csysv_(&uplo, &n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
---|
| 342 | (complex<r_4> *)b.Data(), &ldb,
|
---|
[2554] | 343 | (complex<r_4> *)work, &lwork, &info);
|
---|
| 344 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
| 345 | zsysv_(&uplo, &n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
---|
| 346 | (complex<r_8> *)b.Data(), &ldb,
|
---|
[2572] | 347 | (complex<r_8> *)wkget, &lwork, &info);
|
---|
| 348 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 349 | zsysv_(&uplo, &n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
---|
| 350 | (complex<r_8> *)b.Data(), &ldb,
|
---|
[2554] | 351 | (complex<r_8> *)work, &lwork, &info);
|
---|
| 352 | } else {
|
---|
[2556] | 353 | if(work) delete[] work;
|
---|
[2554] | 354 | delete[] ipiv;
|
---|
| 355 | string tn = typeid(T).name();
|
---|
| 356 | cerr << " LapackServer::LinSolveSym(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 357 | throw TypeMismatchExc("LapackServer::LinSolveSym(a,b) - Unsupported DataType (T)");
|
---|
| 358 | }
|
---|
[2556] | 359 | if(work) delete[] work;
|
---|
[2554] | 360 | delete[] ipiv;
|
---|
| 361 | return(info);
|
---|
| 362 | }
|
---|
| 363 |
|
---|
[2556] | 364 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[1566] | 365 | //! Interface to Lapack least squares solver driver s/d/c/zgels().
|
---|
| 366 | /*! Solves the linear least squares problem defined by an m-by-n matrix
|
---|
[2563] | 367 | \b a and an m element vector \b b , using QR or LQ factorization .
|
---|
[1566] | 368 | A solution \b x to the overdetermined system of linear equations
|
---|
| 369 | b = a * x is computed, minimizing the norm of b-a*x.
|
---|
| 370 | Underdetermined systems (m<n) are not yet handled.
|
---|
| 371 | Inout arrays should have FortranMemory mapping (column packed).
|
---|
| 372 | \param a : input matrix, overwritten on output
|
---|
| 373 | \param b : input-output, input vector b, contains x on exit.
|
---|
| 374 | \return : return code from lapack driver _gels()
|
---|
| 375 | \warning : b is not resized.
|
---|
| 376 | */
|
---|
| 377 | /*
|
---|
| 378 | $CHECK$ - A faire - cas m<n
|
---|
| 379 | If the linear system is underdetermined, the minimum norm
|
---|
| 380 | solution is computed.
|
---|
| 381 | */
|
---|
| 382 |
|
---|
[1494] | 383 | template <class T>
|
---|
| 384 | int LapackServer<T>::LeastSquareSolve(TArray<T>& a, TArray<T> & b)
|
---|
| 385 | {
|
---|
| 386 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
[2561] | 387 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) a Or b NbDimensions() != 2"));
|
---|
[1494] | 388 |
|
---|
| 389 | int_4 rowa = a.RowsKA();
|
---|
| 390 | int_4 cola = a.ColsKA();
|
---|
| 391 | int_4 rowb = b.RowsKA();
|
---|
| 392 | int_4 colb = b.ColsKA();
|
---|
| 393 |
|
---|
| 394 |
|
---|
| 395 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
| 396 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) RowSize(a <> b) "));
|
---|
| 397 |
|
---|
| 398 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
[1566] | 399 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) a Or b Not Column Packed"));
|
---|
[1494] | 400 |
|
---|
[1566] | 401 | if ( a.Size(rowa) < a.Size(cola)) { // $CHECK$ - m<n a changer
|
---|
| 402 | cout << " LapackServer<T>::LeastSquareSolve() - m<n - Not yet implemented for "
|
---|
| 403 | << " underdetermined systems ! " << endl;
|
---|
| 404 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) NRows<NCols - "));
|
---|
| 405 | }
|
---|
[1494] | 406 | int_4 m = a.Size(rowa);
|
---|
| 407 | int_4 n = a.Size(cola);
|
---|
| 408 | int_4 nrhs = b.Size(colb);
|
---|
| 409 |
|
---|
| 410 | int_4 lda = a.Step(cola);
|
---|
| 411 | int_4 ldb = b.Step(colb);
|
---|
| 412 | int_4 info;
|
---|
| 413 |
|
---|
| 414 | int_4 minmn = (m < n) ? m : n;
|
---|
| 415 | int_4 maxmn = (m > n) ? m : n;
|
---|
| 416 | int_4 maxmnrhs = (nrhs > maxmn) ? nrhs : maxmn;
|
---|
| 417 | if (maxmnrhs < 1) maxmnrhs = 1;
|
---|
| 418 |
|
---|
[2572] | 419 | int_4 lwork = -1; //minmn+maxmnrhs*5;
|
---|
| 420 | T * work = NULL;
|
---|
| 421 | T wkget[2];
|
---|
[1494] | 422 |
|
---|
| 423 | char trans = 'N';
|
---|
| 424 |
|
---|
[2572] | 425 | if (typeid(T) == typeid(r_4) ) {
|
---|
[1494] | 426 | sgels_(&trans, &m, &n, &nrhs, (r_4 *)a.Data(), &lda,
|
---|
[2572] | 427 | (r_4 *)b.Data(), &ldb, (r_4 *)wkget, &lwork, &info);
|
---|
| 428 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 429 | sgels_(&trans, &m, &n, &nrhs, (r_4 *)a.Data(), &lda,
|
---|
[1494] | 430 | (r_4 *)b.Data(), &ldb, (r_4 *)work, &lwork, &info);
|
---|
[2572] | 431 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
[1494] | 432 | dgels_(&trans, &m, &n, &nrhs, (r_8 *)a.Data(), &lda,
|
---|
[2572] | 433 | (r_8 *)b.Data(), &ldb, (r_8 *)wkget, &lwork, &info);
|
---|
| 434 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 435 | dgels_(&trans, &m, &n, &nrhs, (r_8 *)a.Data(), &lda,
|
---|
[1494] | 436 | (r_8 *)b.Data(), &ldb, (r_8 *)work, &lwork, &info);
|
---|
[2572] | 437 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
[1494] | 438 | cgels_(&trans, &m, &n, &nrhs, (complex<r_4> *)a.Data(), &lda,
|
---|
[2572] | 439 | (complex<r_4> *)b.Data(), &ldb, (complex<r_4> *)wkget, &lwork, &info);
|
---|
| 440 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 441 | cgels_(&trans, &m, &n, &nrhs, (complex<r_4> *)a.Data(), &lda,
|
---|
[1494] | 442 | (complex<r_4> *)b.Data(), &ldb, (complex<r_4> *)work, &lwork, &info);
|
---|
[2572] | 443 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
[1494] | 444 | zgels_(&trans, &m, &n, &nrhs, (complex<r_8> *)a.Data(), &lda,
|
---|
[2572] | 445 | (complex<r_8> *)b.Data(), &ldb, (complex<r_8> *)wkget, &lwork, &info);
|
---|
| 446 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 447 | zgels_(&trans, &m, &n, &nrhs, (complex<r_8> *)a.Data(), &lda,
|
---|
[1494] | 448 | (complex<r_8> *)b.Data(), &ldb, (complex<r_8> *)work, &lwork, &info);
|
---|
[2572] | 449 | } else {
|
---|
| 450 | if(work) delete [] work; work=NULL;
|
---|
[1494] | 451 | string tn = typeid(T).name();
|
---|
| 452 | cerr << " LapackServer::LeastSquareSolve(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 453 | throw TypeMismatchExc("LapackServer::LeastSquareSolve(a,b) - Unsupported DataType (T)");
|
---|
| 454 | }
|
---|
[2572] | 455 | if(work) delete [] work;
|
---|
[1494] | 456 | return(info);
|
---|
| 457 | }
|
---|
| 458 |
|
---|
[2567] | 459 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
| 460 | //! Interface to Lapack least squares solver driver s/d/c/zgelsd().
|
---|
| 461 | /*! Solves the linear least squares problem defined by an m-by-n matrix
|
---|
| 462 | \b a and an m element vector \b b , using SVD factorization Divide and Conquer.
|
---|
| 463 | Inout arrays should have FortranMemory mapping (column packed).
|
---|
| 464 | \param rcond : definition of zero value (S(i) <= RCOND*S(0) are treated as zero).
|
---|
| 465 | If RCOND < 0, machine precision is used instead.
|
---|
| 466 | \param a : input matrix, overwritten on output
|
---|
| 467 | \param b : input vector b overwritten by solution on output (beware of size changing)
|
---|
| 468 | \param x : output matrix of solutions.
|
---|
[2572] | 469 | \param rank : output the rank of the matrix.
|
---|
[2567] | 470 | \return : return code from lapack driver _gelsd()
|
---|
| 471 | \warning : b is not resized.
|
---|
| 472 | */
|
---|
| 473 | template <class T>
|
---|
| 474 | int LapackServer<T>::LeastSquareSolveSVD_DC(TMatrix<T>& a,TMatrix<T>& b,TVector<r_8>& s,int_4& rank,r_8 rcond)
|
---|
| 475 | {
|
---|
| 476 | if ( ( a.NbDimensions() != 2 ) )
|
---|
| 477 | throw(SzMismatchError("LapackServer::LeastSquareSolveSVD_DC(a,b) a != 2"));
|
---|
| 478 |
|
---|
| 479 | if (!a.IsPacked() || !b.IsPacked())
|
---|
| 480 | throw(SzMismatchError("LapackServer::LeastSquareSolveSVD_DC(a,b) a Or b Not Packed"));
|
---|
[1494] | 481 |
|
---|
[2567] | 482 | int_4 m = a.NRows();
|
---|
| 483 | int_4 n = a.NCols();
|
---|
| 484 |
|
---|
[2572] | 485 | if(b.NRows() != m)
|
---|
[2567] | 486 | throw(SzMismatchError("LapackServer::LeastSquareSolveSVD_DC(a,b) bad matching dim between a and b"));
|
---|
| 487 |
|
---|
| 488 | int_4 nrhs = b.NCols();
|
---|
| 489 | int_4 minmn = (m < n) ? m : n;
|
---|
| 490 | int_4 maxmn = (m > n) ? m : n;
|
---|
| 491 |
|
---|
| 492 | int_4 lda = m;
|
---|
| 493 | int_4 ldb = maxmn;
|
---|
| 494 | int_4 info;
|
---|
| 495 |
|
---|
[2572] | 496 | { // Use {} for automatic des-allocation of "bsave"
|
---|
| 497 | TMatrix<T> bsave(m,nrhs); bsave.SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
[2567] | 498 | bsave = b;
|
---|
| 499 | b.ReSize(maxmn,nrhs); b = (T) 0;
|
---|
[2572] | 500 | for(int i=0;i<m;i++) for(int j=0;j<nrhs;j++) b(i,j) = bsave(i,j);
|
---|
| 501 | } // Use {} for automatic des-allocation of "bsave"
|
---|
[2567] | 502 | s.ReSize(minmn);
|
---|
| 503 |
|
---|
[2572] | 504 | int_4 smlsiz = 25; // Normallement ilaenv_en_C(9,...) renvoie toujours 25
|
---|
[2567] | 505 | if(typeid(T) == typeid(r_4) ) smlsiz = ilaenv_en_C(9,"SGELSD"," ",0,0,0,0);
|
---|
| 506 | else if(typeid(T) == typeid(r_8) ) smlsiz = ilaenv_en_C(9,"DGELSD"," ",0,0,0,0);
|
---|
| 507 | else if(typeid(T) == typeid(complex<r_4>) ) smlsiz = ilaenv_en_C(9,"CGELSD"," ",0,0,0,0);
|
---|
| 508 | else if(typeid(T) == typeid(complex<r_8>) ) smlsiz = ilaenv_en_C(9,"ZGELSD"," ",0,0,0,0);
|
---|
| 509 | if(smlsiz<0) smlsiz = 0;
|
---|
| 510 | r_8 dum = log((r_8)minmn/(r_8)(smlsiz+1.)) / log(2.);
|
---|
| 511 | int_4 nlvl = int_4(dum) + 1; if(nlvl<0) nlvl = 0;
|
---|
| 512 |
|
---|
[2572] | 513 | T * work = NULL;
|
---|
| 514 | int_4 * iwork = NULL;
|
---|
| 515 | int_4 lwork=-1, lrwork;
|
---|
| 516 | T wkget[2];
|
---|
| 517 |
|
---|
[2567] | 518 | if(typeid(T) == typeid(r_4) ) {
|
---|
| 519 | r_4* sloc = new r_4[minmn];
|
---|
| 520 | r_4 srcond = rcond;
|
---|
[2572] | 521 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
---|
[2567] | 522 | sgelsd_(&m,&n,&nrhs,(r_4*)a.Data(),&lda,
|
---|
| 523 | (r_4*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
---|
[2572] | 524 | (r_4*)wkget,&lwork,(int_4*)iwork,&info);
|
---|
| 525 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 526 | sgelsd_(&m,&n,&nrhs,(r_4*)a.Data(),&lda,
|
---|
| 527 | (r_4*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
---|
[2567] | 528 | (r_4*)work,&lwork,(int_4*)iwork,&info);
|
---|
| 529 | for(int_4 i=0;i<minmn;i++) s(i) = sloc[i];
|
---|
[2572] | 530 | delete [] sloc;
|
---|
[2567] | 531 | } else if(typeid(T) == typeid(r_8) ) {
|
---|
[2572] | 532 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
---|
[2567] | 533 | dgelsd_(&m,&n,&nrhs,(r_8*)a.Data(),&lda,
|
---|
| 534 | (r_8*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
---|
[2572] | 535 | (r_8*)wkget,&lwork,(int_4*)iwork,&info);
|
---|
| 536 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 537 | dgelsd_(&m,&n,&nrhs,(r_8*)a.Data(),&lda,
|
---|
| 538 | (r_8*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
---|
[2567] | 539 | (r_8*)work,&lwork,(int_4*)iwork,&info);
|
---|
| 540 | } else if(typeid(T) == typeid(complex<r_4>) ) {
|
---|
[2572] | 541 | // Cf meme remarque que ci-dessous (complex<r_8)
|
---|
| 542 | lrwork = 10*minmn + 2*minmn*smlsiz + 8*minmn*nlvl + 3*smlsiz*nrhs + (smlsiz+1)*(smlsiz+1);
|
---|
| 543 | int_4 lrwork_d = 12*minmn + 2*minmn*smlsiz + 8*minmn*nlvl + minmn*nrhs + (smlsiz+1)*(smlsiz+1);
|
---|
| 544 | if(lrwork_d > lrwork) lrwork = lrwork_d;
|
---|
[2567] | 545 | r_4* rwork = new r_4[lrwork +GARDMEM];
|
---|
[2572] | 546 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
---|
[2567] | 547 | r_4* sloc = new r_4[minmn];
|
---|
| 548 | r_4 srcond = rcond;
|
---|
| 549 | cgelsd_(&m,&n,&nrhs,(complex<r_4>*)a.Data(),&lda,
|
---|
| 550 | (complex<r_4>*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
---|
[2572] | 551 | (complex<r_4>*)wkget,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
---|
| 552 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 553 | cgelsd_(&m,&n,&nrhs,(complex<r_4>*)a.Data(),&lda,
|
---|
| 554 | (complex<r_4>*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
---|
[2567] | 555 | (complex<r_4>*)work,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
---|
| 556 | for(int_4 i=0;i<minmn;i++) s(i) = sloc[i];
|
---|
[2572] | 557 | delete [] sloc; delete [] rwork;
|
---|
[2567] | 558 | } else if(typeid(T) == typeid(complex<r_8>) ) {
|
---|
[2572] | 559 | // CMV: Bizarrement, la formule donnee dans zgelsd() plante pour des N grands (500)
|
---|
| 560 | // On prend (par analogie) la formule pour "lwork" de dgelsd()
|
---|
| 561 | lrwork = 10*minmn + 2*minmn*smlsiz + 8*minmn*nlvl + 3*smlsiz*nrhs + (smlsiz+1)*(smlsiz+1);
|
---|
| 562 | int_4 lrwork_d = 12*minmn + 2*minmn*smlsiz + 8*minmn*nlvl + minmn*nrhs + (smlsiz+1)*(smlsiz+1);
|
---|
| 563 | if(lrwork_d > lrwork) lrwork = lrwork_d;
|
---|
[2567] | 564 | r_8* rwork = new r_8[lrwork +GARDMEM];
|
---|
[2572] | 565 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
---|
[2567] | 566 | zgelsd_(&m,&n,&nrhs,(complex<r_8>*)a.Data(),&lda,
|
---|
[2572] | 567 | (complex<r_8>*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
---|
| 568 | (complex<r_8>*)wkget,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
---|
| 569 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 570 | zgelsd_(&m,&n,&nrhs,(complex<r_8>*)a.Data(),&lda,
|
---|
| 571 | (complex<r_8>*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
---|
[2567] | 572 | (complex<r_8>*)work,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
---|
[2572] | 573 | delete [] rwork;
|
---|
[2567] | 574 | } else {
|
---|
[2572] | 575 | if(work) delete [] work; work=NULL;
|
---|
| 576 | if(iwork) delete [] iwork; iwork=NULL;
|
---|
[2567] | 577 | string tn = typeid(T).name();
|
---|
| 578 | cerr << " LapackServer::LeastSquareSolveSVD_DC(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 579 | throw TypeMismatchExc("LapackServer::LeastSquareSolveSVD_DC(a,b) - Unsupported DataType (T)");
|
---|
| 580 | }
|
---|
| 581 |
|
---|
[2572] | 582 | if(work) delete [] work; if(iwork) delete [] iwork;
|
---|
[2567] | 583 | return(info);
|
---|
| 584 | }
|
---|
| 585 |
|
---|
| 586 |
|
---|
[2556] | 587 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[1424] | 588 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
| 589 | /*! Computes the vector of singular values of \b a. Input arrays
|
---|
| 590 | should have FortranMemoryMapping (column packed).
|
---|
| 591 | \param a : input m-by-n matrix
|
---|
| 592 | \param s : Vector of min(m,n) singular values (descending order)
|
---|
| 593 | \return : return code from lapack driver _gesvd()
|
---|
| 594 | */
|
---|
| 595 |
|
---|
[1342] | 596 | template <class T>
|
---|
| 597 | int LapackServer<T>::SVD(TArray<T>& a, TArray<T> & s)
|
---|
| 598 | {
|
---|
| 599 | return (SVDDriver(a, s, NULL, NULL) );
|
---|
| 600 | }
|
---|
| 601 |
|
---|
[1424] | 602 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
| 603 | /*! Computes the vector of singular values of \b a, as well as
|
---|
| 604 | right and left singular vectors of \b a.
|
---|
| 605 | \f[
|
---|
| 606 | A = U \Sigma V^T , ( A = U \Sigma V^H \ complex)
|
---|
| 607 | \f]
|
---|
| 608 | \f[
|
---|
| 609 | A v_i = \sigma_i u_i \ and A^T u_i = \sigma_i v_i \ (A^H \ complex)
|
---|
| 610 | \f]
|
---|
| 611 | U and V are orthogonal (unitary) matrices.
|
---|
[2572] | 612 | \param a : input m-by-n matrix (in FortranMemoryMapping)
|
---|
[1424] | 613 | \param s : Vector of min(m,n) singular values (descending order)
|
---|
[2572] | 614 | \param u : m-by-m Matrix of left singular vectors
|
---|
| 615 | \param vt : Transpose of right singular vectors (n-by-n matrix).
|
---|
[1424] | 616 | \return : return code from lapack driver _gesvd()
|
---|
| 617 | */
|
---|
[1342] | 618 | template <class T>
|
---|
| 619 | int LapackServer<T>::SVD(TArray<T>& a, TArray<T> & s, TArray<T> & u, TArray<T> & vt)
|
---|
| 620 | {
|
---|
| 621 | return (SVDDriver(a, s, &u, &vt) );
|
---|
| 622 | }
|
---|
| 623 |
|
---|
[1424] | 624 |
|
---|
| 625 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
[1342] | 626 | template <class T>
|
---|
| 627 | int LapackServer<T>::SVDDriver(TArray<T>& a, TArray<T> & s, TArray<T>* up, TArray<T>* vtp)
|
---|
| 628 | {
|
---|
| 629 | if ( ( a.NbDimensions() != 2 ) )
|
---|
[2561] | 630 | throw(SzMismatchError("LapackServer::SVDDriver(a, ...) a.NbDimensions() != 2"));
|
---|
[1342] | 631 |
|
---|
| 632 | int_4 rowa = a.RowsKA();
|
---|
| 633 | int_4 cola = a.ColsKA();
|
---|
| 634 |
|
---|
| 635 | if ( !a.IsPacked(rowa) )
|
---|
[2561] | 636 | throw(SzMismatchError("LapackServer::SVDDriver(a, ...) a Not Column Packed "));
|
---|
[1342] | 637 |
|
---|
| 638 | int_4 m = a.Size(rowa);
|
---|
| 639 | int_4 n = a.Size(cola);
|
---|
| 640 | int_4 maxmn = (m > n) ? m : n;
|
---|
| 641 | int_4 minmn = (m < n) ? m : n;
|
---|
| 642 |
|
---|
| 643 | char jobu, jobvt;
|
---|
| 644 | jobu = 'N';
|
---|
| 645 | jobvt = 'N';
|
---|
| 646 |
|
---|
| 647 | sa_size_t sz[2];
|
---|
| 648 | if ( up != NULL) {
|
---|
| 649 | if ( dynamic_cast< TVector<T> * > (vtp) )
|
---|
[2561] | 650 | throw( TypeMismatchExc("LapackServer::SVDDriver() Wrong type (=TVector<T>) for u !") );
|
---|
[1342] | 651 | up->SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
| 652 | sz[0] = sz[1] = m;
|
---|
| 653 | up->ReSize(2, sz );
|
---|
| 654 | jobu = 'A';
|
---|
| 655 | }
|
---|
| 656 | else {
|
---|
| 657 | up = new TMatrix<T>(1,1);
|
---|
| 658 | jobu = 'N';
|
---|
| 659 | }
|
---|
| 660 | if ( vtp != NULL) {
|
---|
| 661 | if ( dynamic_cast< TVector<T> * > (vtp) )
|
---|
[2561] | 662 | throw( TypeMismatchExc("LapackServer::SVDDriver() Wrong type (=TVector<T>) for vt !") );
|
---|
[1342] | 663 | vtp->SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
| 664 | sz[0] = sz[1] = n;
|
---|
| 665 | vtp->ReSize(2, sz );
|
---|
| 666 | jobvt = 'A';
|
---|
| 667 | }
|
---|
| 668 | else {
|
---|
| 669 | vtp = new TMatrix<T>(1,1);
|
---|
| 670 | jobvt = 'N';
|
---|
| 671 | }
|
---|
| 672 |
|
---|
| 673 | TVector<T> *vs = dynamic_cast< TVector<T> * > (&s);
|
---|
| 674 | if (vs) vs->ReSize(minmn);
|
---|
| 675 | else {
|
---|
| 676 | TMatrix<T> *ms = dynamic_cast< TMatrix<T> * > (&s);
|
---|
| 677 | if (ms) ms->ReSize(minmn,1);
|
---|
| 678 | else {
|
---|
| 679 | sz[0] = minmn; sz[1] = 1;
|
---|
| 680 | s.ReSize(1, sz);
|
---|
| 681 | }
|
---|
| 682 | }
|
---|
| 683 |
|
---|
| 684 | int_4 lda = a.Step(a.ColsKA());
|
---|
| 685 | int_4 ldu = up->Step(up->ColsKA());
|
---|
| 686 | int_4 ldvt = vtp->Step(vtp->ColsKA());
|
---|
[2567] | 687 | int_4 info;
|
---|
[1342] | 688 |
|
---|
[2572] | 689 | int_4 lwork = -1; // maxmn*5 *wspace_size_factor;
|
---|
| 690 | T * work = NULL; // = new T[lwork];
|
---|
| 691 | T wkget[2];
|
---|
[1342] | 692 |
|
---|
[2559] | 693 | if (typeid(T) == typeid(r_4) ) {
|
---|
[1342] | 694 | sgesvd_(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda,
|
---|
| 695 | (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt,
|
---|
[2572] | 696 | (r_4 *)wkget, &lwork, &info);
|
---|
| 697 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 698 | sgesvd_(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda,
|
---|
| 699 | (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt,
|
---|
[1342] | 700 | (r_4 *)work, &lwork, &info);
|
---|
[2559] | 701 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
[1342] | 702 | dgesvd_(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda,
|
---|
| 703 | (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt,
|
---|
[2572] | 704 | (r_8 *)wkget, &lwork, &info);
|
---|
| 705 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 706 | dgesvd_(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda,
|
---|
| 707 | (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt,
|
---|
[1342] | 708 | (r_8 *)work, &lwork, &info);
|
---|
[2559] | 709 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
[2567] | 710 | r_4 * rwork = new r_4[5*minmn +GARDMEM];
|
---|
[2559] | 711 | r_4 * sloc = new r_4[minmn];
|
---|
[1342] | 712 | cgesvd_(&jobu, &jobvt, &m, &n, (complex<r_4> *)a.Data(), &lda,
|
---|
[2559] | 713 | (r_4 *)sloc, (complex<r_4> *) up->Data(), &ldu,
|
---|
[1342] | 714 | (complex<r_4> *)vtp->Data(), &ldvt,
|
---|
[2572] | 715 | (complex<r_4> *)wkget, &lwork, (r_4 *)rwork, &info);
|
---|
| 716 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 717 | cgesvd_(&jobu, &jobvt, &m, &n, (complex<r_4> *)a.Data(), &lda,
|
---|
| 718 | (r_4 *)sloc, (complex<r_4> *) up->Data(), &ldu,
|
---|
| 719 | (complex<r_4> *)vtp->Data(), &ldvt,
|
---|
[2559] | 720 | (complex<r_4> *)work, &lwork, (r_4 *)rwork, &info);
|
---|
| 721 | for(int_4 i=0;i<minmn;i++) s[i] = sloc[i];
|
---|
| 722 | delete [] rwork; delete [] sloc;
|
---|
| 723 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
[2567] | 724 | r_8 * rwork = new r_8[5*minmn +GARDMEM];
|
---|
[2559] | 725 | r_8 * sloc = new r_8[minmn];
|
---|
[1342] | 726 | zgesvd_(&jobu, &jobvt, &m, &n, (complex<r_8> *)a.Data(), &lda,
|
---|
[2559] | 727 | (r_8 *)sloc, (complex<r_8> *) up->Data(), &ldu,
|
---|
[1342] | 728 | (complex<r_8> *)vtp->Data(), &ldvt,
|
---|
[2572] | 729 | (complex<r_8> *)wkget, &lwork, (r_8 *)rwork, &info);
|
---|
| 730 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 731 | zgesvd_(&jobu, &jobvt, &m, &n, (complex<r_8> *)a.Data(), &lda,
|
---|
| 732 | (r_8 *)sloc, (complex<r_8> *) up->Data(), &ldu,
|
---|
| 733 | (complex<r_8> *)vtp->Data(), &ldvt,
|
---|
[2559] | 734 | (complex<r_8> *)work, &lwork, (r_8 *)rwork, &info);
|
---|
| 735 | for(int_4 i=0;i<minmn;i++) s[i] = sloc[i];
|
---|
| 736 | delete [] rwork; delete [] sloc;
|
---|
| 737 | } else {
|
---|
[2572] | 738 | if(work) delete [] work; work=NULL;
|
---|
[1342] | 739 | if (jobu == 'N') delete up;
|
---|
| 740 | if (jobvt == 'N') delete vtp;
|
---|
| 741 | string tn = typeid(T).name();
|
---|
| 742 | cerr << " LapackServer::SVDDriver(...) - Unsupported DataType T = " << tn << endl;
|
---|
[2561] | 743 | throw TypeMismatchExc("LapackServer::SVDDriver(a,b) - Unsupported DataType (T)");
|
---|
[1342] | 744 | }
|
---|
| 745 |
|
---|
[2572] | 746 | if(work) delete [] work;
|
---|
[1342] | 747 | if (jobu == 'N') delete up;
|
---|
| 748 | if (jobvt == 'N') delete vtp;
|
---|
| 749 | return(info);
|
---|
| 750 | }
|
---|
| 751 |
|
---|
[2556] | 752 |
|
---|
[2561] | 753 | //! Interface to Lapack SVD driver s/d/c/zgesdd().
|
---|
| 754 | /*! Same as SVD but with Divide and Conquer method */
|
---|
| 755 | template <class T>
|
---|
[2563] | 756 | int LapackServer<T>::SVD_DC(TMatrix<T>& a, TVector<r_8>& s, TMatrix<T>& u, TMatrix<T>& vt)
|
---|
[2561] | 757 | {
|
---|
| 758 |
|
---|
| 759 | if ( !a.IsPacked() )
|
---|
| 760 | throw(SzMismatchError("LapackServer::SVD_DC(a, ...) a Not Packed "));
|
---|
| 761 |
|
---|
| 762 | int_4 m = a.NRows();
|
---|
| 763 | int_4 n = a.NCols();
|
---|
| 764 | int_4 maxmn = (m > n) ? m : n;
|
---|
| 765 | int_4 minmn = (m < n) ? m : n;
|
---|
| 766 | int_4 supermax = 4*minmn*minmn+4*minmn; if(maxmn>supermax) supermax=maxmn;
|
---|
| 767 |
|
---|
| 768 | char jobz = 'A';
|
---|
| 769 |
|
---|
| 770 | s.ReSize(minmn);
|
---|
| 771 | u.ReSize(m,m);
|
---|
| 772 | vt.ReSize(n,n);
|
---|
| 773 |
|
---|
[2572] | 774 | int_4 lda = m;
|
---|
| 775 | int_4 ldu = m;
|
---|
| 776 | int_4 ldvt = n;
|
---|
| 777 | int_4 info;
|
---|
| 778 | int_4 lwork=-1;
|
---|
| 779 | T * work = NULL;
|
---|
| 780 | int_4 * iwork = NULL;
|
---|
| 781 | T wkget[2];
|
---|
| 782 |
|
---|
[2561] | 783 | if(typeid(T) == typeid(r_4) ) {
|
---|
[2563] | 784 | r_4* sloc = new r_4[minmn];
|
---|
[2572] | 785 | iwork = new int_4[8*minmn +GARDMEM];
|
---|
[2561] | 786 | sgesdd_(&jobz,&m,&n,(r_4*)a.Data(),&lda,
|
---|
[2563] | 787 | (r_4*)sloc,(r_4*)u.Data(),&ldu,(r_4*)vt.Data(),&ldvt,
|
---|
[2572] | 788 | (r_4*)wkget,&lwork,(int_4*)iwork,&info);
|
---|
| 789 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 790 | sgesdd_(&jobz,&m,&n,(r_4*)a.Data(),&lda,
|
---|
| 791 | (r_4*)sloc,(r_4*)u.Data(),&ldu,(r_4*)vt.Data(),&ldvt,
|
---|
[2561] | 792 | (r_4*)work,&lwork,(int_4*)iwork,&info);
|
---|
[2563] | 793 | for(int_4 i=0;i<minmn;i++) s(i) = (r_8) sloc[i];
|
---|
[2572] | 794 | delete [] sloc;
|
---|
[2561] | 795 | } else if(typeid(T) == typeid(r_8) ) {
|
---|
[2572] | 796 | iwork = new int_4[8*minmn +GARDMEM];
|
---|
[2561] | 797 | dgesdd_(&jobz,&m,&n,(r_8*)a.Data(),&lda,
|
---|
| 798 | (r_8*)s.Data(),(r_8*)u.Data(),&ldu,(r_8*)vt.Data(),&ldvt,
|
---|
[2572] | 799 | (r_8*)wkget,&lwork,(int_4*)iwork,&info);
|
---|
| 800 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 801 | dgesdd_(&jobz,&m,&n,(r_8*)a.Data(),&lda,
|
---|
| 802 | (r_8*)s.Data(),(r_8*)u.Data(),&ldu,(r_8*)vt.Data(),&ldvt,
|
---|
[2561] | 803 | (r_8*)work,&lwork,(int_4*)iwork,&info);
|
---|
| 804 | } else if(typeid(T) == typeid(complex<r_4>) ) {
|
---|
| 805 | r_4* sloc = new r_4[minmn];
|
---|
[2567] | 806 | r_4* rwork = new r_4[5*minmn*minmn+5*minmn +GARDMEM];
|
---|
[2572] | 807 | iwork = new int_4[8*minmn +GARDMEM];
|
---|
[2561] | 808 | cgesdd_(&jobz,&m,&n,(complex<r_4>*)a.Data(),&lda,
|
---|
| 809 | (r_4*)sloc,(complex<r_4>*)u.Data(),&ldu,(complex<r_4>*)vt.Data(),&ldvt,
|
---|
[2572] | 810 | (complex<r_4>*)wkget,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
---|
| 811 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
---|
| 812 | cgesdd_(&jobz,&m,&n,(complex<r_4>*)a.Data(),&lda,
|
---|
| 813 | (r_4*)sloc,(complex<r_4>*)u.Data(),&ldu,(complex<r_4>*)vt.Data(),&ldvt,
|
---|
[2561] | 814 | (complex<r_4>*)work,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
---|
[2563] | 815 | for(int_4 i=0;i<minmn;i++) s(i) = (r_8) sloc[i];
|
---|
[2572] | 816 | delete [] sloc; delete [] rwork;
|
---|
[2561] | 817 | } else if(typeid(T) == typeid(complex<r_8>) ) {
|
---|
[2567] | 818 | r_8* rwork = new r_8[5*minmn*minmn+5*minmn +GARDMEM];
|
---|
[2572] | 819 | iwork = new int_4[8*minmn +GARDMEM];
|
---|
[2561] | 820 | zgesdd_(&jobz,&m,&n,(complex<r_8>*)a.Data(),&lda,
|
---|
[2563] | 821 | (r_8*)s.Data(),(complex<r_8>*)u.Data(),&ldu,(complex<r_8>*)vt.Data(),&ldvt,
|
---|
[2572] | 822 | (complex<r_8>*)wkget,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
---|
| 823 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
---|
| 824 | zgesdd_(&jobz,&m,&n,(complex<r_8>*)a.Data(),&lda,
|
---|
| 825 | (r_8*)s.Data(),(complex<r_8>*)u.Data(),&ldu,(complex<r_8>*)vt.Data(),&ldvt,
|
---|
[2561] | 826 | (complex<r_8>*)work,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
---|
[2572] | 827 | delete [] rwork;
|
---|
[2561] | 828 | } else {
|
---|
[2572] | 829 | if(work) delete [] work; work=NULL;
|
---|
| 830 | if(iwork) delete [] iwork; iwork=NULL;
|
---|
[2561] | 831 | string tn = typeid(T).name();
|
---|
| 832 | cerr << " LapackServer::SVD_DC(...) - Unsupported DataType T = " << tn << endl;
|
---|
| 833 | throw TypeMismatchExc("LapackServer::SVD_DC - Unsupported DataType (T)");
|
---|
| 834 | }
|
---|
| 835 |
|
---|
[2572] | 836 | if(work) delete [] work; if(iwork) delete [] iwork;
|
---|
[2561] | 837 | return(info);
|
---|
| 838 | }
|
---|
| 839 |
|
---|
| 840 |
|
---|
[2556] | 841 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
| 842 | /*! Computes the eigen values and eigen vectors of a symetric (or hermitian) matrix \b a.
|
---|
| 843 | Input arrays should have FortranMemoryMapping (column packed).
|
---|
| 844 | \param a : input symetric (or hermitian) n-by-n matrix
|
---|
| 845 | \param b : Vector of eigenvalues (descending order)
|
---|
| 846 | \param eigenvector : if true compute eigenvectors, if not only eigenvalues
|
---|
| 847 | \param a : on return array of eigenvectors (same order than eval, one vector = one column)
|
---|
[2561] | 848 | \return : return code from lapack driver
|
---|
[2556] | 849 | */
|
---|
| 850 |
|
---|
| 851 | template <class T>
|
---|
| 852 | int LapackServer<T>::LapackEigenSym(TArray<T>& a, TVector<r_8>& b, bool eigenvector)
|
---|
| 853 | {
|
---|
| 854 | if ( a.NbDimensions() != 2 )
|
---|
| 855 | throw(SzMismatchError("LapackServer::LapackEigenSym(a,b) a NbDimensions() != 2"));
|
---|
| 856 | int_4 rowa = a.RowsKA();
|
---|
| 857 | int_4 cola = a.ColsKA();
|
---|
| 858 | if ( a.Size(rowa) != a.Size(cola))
|
---|
| 859 | throw(SzMismatchError("LapackServer::LapackEigenSym(a,b) a Not a square Array"));
|
---|
| 860 | if (!a.IsPacked(rowa))
|
---|
| 861 | throw(SzMismatchError("LapackServer::LapackEigenSym(a,b) a Not Column Packed"));
|
---|
| 862 |
|
---|
[2561] | 863 | char uplo='U';
|
---|
[2556] | 864 | char jobz='N'; if(eigenvector) jobz='V';
|
---|
| 865 |
|
---|
| 866 | int_4 n = a.Size(rowa);
|
---|
| 867 | int_4 lda = a.Step(cola);
|
---|
| 868 | int_4 info = 0;
|
---|
[2572] | 869 | int_4 lwork = -1;
|
---|
| 870 | T * work = NULL;
|
---|
| 871 | T wkget[2];
|
---|
[2556] | 872 |
|
---|
| 873 | b.ReSize(n); b = 0.;
|
---|
| 874 |
|
---|
| 875 | if (typeid(T) == typeid(r_4) ) {
|
---|
| 876 | r_4* w = new r_4[n];
|
---|
[2572] | 877 | ssyev_(&jobz,&uplo,&n,(r_4 *)a.Data(),&lda,(r_4 *)w,(r_4 *)wkget,&lwork,&info);
|
---|
| 878 | lwork = type2i4(&wkget[0],4); /* 3*n-1;*/ work = new T[lwork +GARDMEM];
|
---|
[2561] | 879 | ssyev_(&jobz,&uplo,&n,(r_4 *)a.Data(),&lda,(r_4 *)w,(r_4 *)work,&lwork,&info);
|
---|
[2556] | 880 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
---|
[2572] | 881 | delete [] w;
|
---|
[2556] | 882 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
| 883 | r_8* w = new r_8[n];
|
---|
[2572] | 884 | dsyev_(&jobz,&uplo,&n,(r_8 *)a.Data(),&lda,(r_8 *)w,(r_8 *)wkget,&lwork,&info);
|
---|
| 885 | lwork = type2i4(&wkget[0],8); /* 3*n-1;*/ work = new T[lwork +GARDMEM];
|
---|
[2561] | 886 | dsyev_(&jobz,&uplo,&n,(r_8 *)a.Data(),&lda,(r_8 *)w,(r_8 *)work,&lwork,&info);
|
---|
[2556] | 887 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
---|
[2572] | 888 | delete [] w;
|
---|
[2556] | 889 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
[2567] | 890 | r_4* rwork = new r_4[3*n-2 +GARDMEM]; r_4* w = new r_4[n];
|
---|
[2561] | 891 | cheev_(&jobz,&uplo,&n,(complex<r_4> *)a.Data(),&lda,(r_4 *)w
|
---|
[2572] | 892 | ,(complex<r_4> *)wkget,&lwork,(r_4 *)rwork,&info);
|
---|
| 893 | lwork = type2i4(&wkget[0],4); /* 2*n-1;*/ work = new T[lwork +GARDMEM];
|
---|
| 894 | cheev_(&jobz,&uplo,&n,(complex<r_4> *)a.Data(),&lda,(r_4 *)w
|
---|
[2556] | 895 | ,(complex<r_4> *)work,&lwork,(r_4 *)rwork,&info);
|
---|
| 896 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
---|
[2572] | 897 | delete [] rwork; delete [] w;
|
---|
[2556] | 898 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
[2567] | 899 | r_8* rwork = new r_8[3*n-2 +GARDMEM]; r_8* w = new r_8[n];
|
---|
[2561] | 900 | zheev_(&jobz,&uplo,&n,(complex<r_8> *)a.Data(),&lda,(r_8 *)w
|
---|
[2572] | 901 | ,(complex<r_8> *)wkget,&lwork,(r_8 *)rwork,&info);
|
---|
| 902 | lwork = type2i4(&wkget[0],8); /* 2*n-1;*/ work = new T[lwork +GARDMEM];
|
---|
| 903 | zheev_(&jobz,&uplo,&n,(complex<r_8> *)a.Data(),&lda,(r_8 *)w
|
---|
[2556] | 904 | ,(complex<r_8> *)work,&lwork,(r_8 *)rwork,&info);
|
---|
| 905 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
---|
[2572] | 906 | delete [] rwork; delete [] w;
|
---|
[2556] | 907 | } else {
|
---|
[2572] | 908 | if(work) delete [] work; work=NULL;
|
---|
[2556] | 909 | string tn = typeid(T).name();
|
---|
| 910 | cerr << " LapackServer::LapackEigenSym(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 911 | throw TypeMismatchExc("LapackServer::LapackEigenSym(a,b) - Unsupported DataType (T)");
|
---|
| 912 | }
|
---|
| 913 |
|
---|
[2572] | 914 | if(work) delete [] work;
|
---|
[2556] | 915 | return(info);
|
---|
| 916 | }
|
---|
| 917 |
|
---|
| 918 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
| 919 | /*! Computes the eigen values and eigen vectors of a general squared matrix \b a.
|
---|
| 920 | Input arrays should have FortranMemoryMapping (column packed).
|
---|
| 921 | \param a : input general n-by-n matrix
|
---|
| 922 | \param eval : Vector of eigenvalues (complex double precision)
|
---|
| 923 | \param evec : Matrix of eigenvector (same order than eval, one vector = one column)
|
---|
| 924 | \param eigenvector : if true compute (right) eigenvectors, if not only eigenvalues
|
---|
| 925 | \param a : on return array of eigenvectors
|
---|
[2561] | 926 | \return : return code from lapack driver
|
---|
[2556] | 927 | \verbatim
|
---|
| 928 | eval : contains the computed eigenvalues.
|
---|
| 929 | --- For real matrices "a" :
|
---|
| 930 | Complex conjugate pairs of eigenvalues appear consecutively
|
---|
| 931 | with the eigenvalue having the positive imaginary part first.
|
---|
| 932 | evec : the right eigenvectors v(j) are stored one after another
|
---|
| 933 | in the columns of evec, in the same order as their eigenvalues.
|
---|
| 934 | --- For real matrices "a" :
|
---|
| 935 | If the j-th eigenvalue is real, then v(j) = evec(:,j),
|
---|
| 936 | the j-th column of evec.
|
---|
| 937 | If the j-th and (j+1)-st eigenvalues form a complex
|
---|
| 938 | conjugate pair, then v(j) = evec(:,j) + i*evec(:,j+1) and
|
---|
| 939 | v(j+1) = evec(:,j) - i*evec(:,j+1).
|
---|
| 940 | \endverbatim
|
---|
| 941 | */
|
---|
| 942 |
|
---|
| 943 | template <class T>
|
---|
| 944 | int LapackServer<T>::LapackEigen(TArray<T>& a, TVector< complex<r_8> >& eval, TMatrix<T>& evec, bool eigenvector)
|
---|
| 945 | {
|
---|
| 946 | if ( a.NbDimensions() != 2 )
|
---|
| 947 | throw(SzMismatchError("LapackServer::LapackEigen(a,b) a NbDimensions() != 2"));
|
---|
| 948 | int_4 rowa = a.RowsKA();
|
---|
| 949 | int_4 cola = a.ColsKA();
|
---|
| 950 | if ( a.Size(rowa) != a.Size(cola))
|
---|
| 951 | throw(SzMismatchError("LapackServer::LapackEigen(a,b) a Not a square Array"));
|
---|
| 952 | if (!a.IsPacked(rowa))
|
---|
| 953 | throw(SzMismatchError("LapackServer::LapackEigen(a,b) a Not Column Packed"));
|
---|
| 954 |
|
---|
[2561] | 955 | char jobvl = 'N';
|
---|
[2556] | 956 | char jobvr = 'N'; if(eigenvector) jobvr='V';
|
---|
| 957 |
|
---|
| 958 | int_4 n = a.Size(rowa);
|
---|
| 959 | int_4 lda = a.Step(cola);
|
---|
| 960 | int_4 info = 0;
|
---|
| 961 |
|
---|
| 962 | eval.ReSize(n); eval = complex<r_8>(0.,0.);
|
---|
| 963 | if(eigenvector) {evec.ReSize(n,n); evec = (T) 0.;}
|
---|
| 964 | int_4 ldvr = n, ldvl = 1;
|
---|
| 965 |
|
---|
[2572] | 966 | int_4 lwork = -1;
|
---|
| 967 | T * work = NULL;
|
---|
| 968 | T wkget[2];
|
---|
| 969 |
|
---|
[2556] | 970 | if (typeid(T) == typeid(r_4) ) {
|
---|
| 971 | r_4* wr = new r_4[n]; r_4* wi = new r_4[n]; r_4* vl = NULL;
|
---|
[2561] | 972 | sgeev_(&jobvl,&jobvr,&n,(r_4 *)a.Data(),&lda,(r_4 *)wr,(r_4 *)wi,
|
---|
[2556] | 973 | (r_4 *)vl,&ldvl,(r_4 *)evec.Data(),&ldvr,
|
---|
[2572] | 974 | (r_4 *)wkget,&lwork,&info);
|
---|
| 975 | lwork = type2i4(&wkget[0],4); /* 4*n;*/ work = new T[lwork +GARDMEM];
|
---|
| 976 | sgeev_(&jobvl,&jobvr,&n,(r_4 *)a.Data(),&lda,(r_4 *)wr,(r_4 *)wi,
|
---|
| 977 | (r_4 *)vl,&ldvl,(r_4 *)evec.Data(),&ldvr,
|
---|
[2556] | 978 | (r_4 *)work,&lwork,&info);
|
---|
| 979 | if(info==0) for(int i=0;i<n;i++) eval(i) = complex<r_8>(wr[i],wi[i]);
|
---|
[2572] | 980 | delete [] wr; delete [] wi;
|
---|
[2556] | 981 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
| 982 | r_8* wr = new r_8[n]; r_8* wi = new r_8[n]; r_8* vl = NULL;
|
---|
[2561] | 983 | dgeev_(&jobvl,&jobvr,&n,(r_8 *)a.Data(),&lda,(r_8 *)wr,(r_8 *)wi,
|
---|
[2556] | 984 | (r_8 *)vl,&ldvl,(r_8 *)evec.Data(),&ldvr,
|
---|
[2572] | 985 | (r_8 *)wkget,&lwork,&info);
|
---|
| 986 | lwork = type2i4(&wkget[0],8); /* 4*n;*/ work = new T[lwork +GARDMEM];
|
---|
| 987 | dgeev_(&jobvl,&jobvr,&n,(r_8 *)a.Data(),&lda,(r_8 *)wr,(r_8 *)wi,
|
---|
| 988 | (r_8 *)vl,&ldvl,(r_8 *)evec.Data(),&ldvr,
|
---|
[2556] | 989 | (r_8 *)work,&lwork,&info);
|
---|
| 990 | if(info==0) for(int i=0;i<n;i++) eval(i) = complex<r_8>(wr[i],wi[i]);
|
---|
[2572] | 991 | delete [] wr; delete [] wi;
|
---|
[2556] | 992 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
[2567] | 993 | r_4* rwork = new r_4[2*n +GARDMEM]; r_4* vl = NULL; TVector< complex<r_4> > w(n);
|
---|
[2561] | 994 | cgeev_(&jobvl,&jobvr,&n,(complex<r_4> *)a.Data(),&lda,(complex<r_4> *)w.Data(),
|
---|
[2556] | 995 | (complex<r_4> *)vl,&ldvl,(complex<r_4> *)evec.Data(),&ldvr,
|
---|
[2572] | 996 | (complex<r_4> *)wkget,&lwork,(r_4 *)rwork,&info);
|
---|
| 997 | lwork = type2i4(&wkget[0],4); /* 2*n;*/ work = new T[lwork +GARDMEM];
|
---|
| 998 | cgeev_(&jobvl,&jobvr,&n,(complex<r_4> *)a.Data(),&lda,(complex<r_4> *)w.Data(),
|
---|
| 999 | (complex<r_4> *)vl,&ldvl,(complex<r_4> *)evec.Data(),&ldvr,
|
---|
[2556] | 1000 | (complex<r_4> *)work,&lwork,(r_4 *)rwork,&info);
|
---|
| 1001 | if(info==0) for(int i=0;i<n;i++) eval(i) = w(i);
|
---|
[2572] | 1002 | delete [] rwork;
|
---|
[2556] | 1003 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
[2567] | 1004 | r_8* rwork = new r_8[2*n +GARDMEM]; r_8* vl = NULL;
|
---|
[2561] | 1005 | zgeev_(&jobvl,&jobvr,&n,(complex<r_8> *)a.Data(),&lda,(complex<r_8> *)eval.Data(),
|
---|
[2556] | 1006 | (complex<r_8> *)vl,&ldvl,(complex<r_8> *)evec.Data(),&ldvr,
|
---|
[2572] | 1007 | (complex<r_8> *)wkget,&lwork,(r_8 *)rwork,&info);
|
---|
| 1008 | lwork = type2i4(&wkget[0],8); /* 2*n;*/ work = new T[lwork +GARDMEM];
|
---|
| 1009 | zgeev_(&jobvl,&jobvr,&n,(complex<r_8> *)a.Data(),&lda,(complex<r_8> *)eval.Data(),
|
---|
| 1010 | (complex<r_8> *)vl,&ldvl,(complex<r_8> *)evec.Data(),&ldvr,
|
---|
[2556] | 1011 | (complex<r_8> *)work,&lwork,(r_8 *)rwork,&info);
|
---|
[2572] | 1012 | delete [] rwork;
|
---|
[2556] | 1013 | } else {
|
---|
[2572] | 1014 | if(work) delete [] work; work=NULL;
|
---|
[2556] | 1015 | string tn = typeid(T).name();
|
---|
| 1016 | cerr << " LapackServer::LapackEigen(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
| 1017 | throw TypeMismatchExc("LapackServer::LapackEigen(a,b) - Unsupported DataType (T)");
|
---|
| 1018 | }
|
---|
| 1019 |
|
---|
[2572] | 1020 | if(work) delete [] work;
|
---|
[2556] | 1021 | return(info);
|
---|
| 1022 | }
|
---|
| 1023 |
|
---|
| 1024 |
|
---|
| 1025 |
|
---|
| 1026 |
|
---|
| 1027 |
|
---|
| 1028 | ////////////////////////////////////////////////////////////////////////////////////
|
---|
[775] | 1029 | void rztest_lapack(TArray<r_4>& aa, TArray<r_4>& bb)
|
---|
| 1030 | {
|
---|
| 1031 | if ( aa.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a Matrix"));
|
---|
| 1032 | if ( aa.SizeX() != aa.SizeY()) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a square Matrix"));
|
---|
| 1033 | if ( bb.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a Matrix"));
|
---|
[788] | 1034 | if ( bb.SizeX() != aa.SizeX() ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A <> B "));
|
---|
[775] | 1035 | if ( !bb.IsPacked() || !bb.IsPacked() )
|
---|
| 1036 | throw(SzMismatchError("rztest_lapack(TMatrix<r_4> Not packed A or B "));
|
---|
| 1037 |
|
---|
[788] | 1038 | int_4 n = aa.SizeX();
|
---|
| 1039 | int_4 nrhs = bb.SizeY();
|
---|
[775] | 1040 | int_4 lda = n;
|
---|
[788] | 1041 | int_4 ldb = bb.SizeX();
|
---|
[775] | 1042 | int_4 info;
|
---|
| 1043 | int_4* ipiv = new int_4[n];
|
---|
| 1044 | sgesv_(&n, &nrhs, aa.Data(), &lda, ipiv, bb.Data(), &ldb, &info);
|
---|
[814] | 1045 | delete[] ipiv;
|
---|
[775] | 1046 | cout << "rztest_lapack/Info= " << info << endl;
|
---|
| 1047 | cout << aa << "\n" << bb << endl;
|
---|
| 1048 | return;
|
---|
| 1049 | }
|
---|
[814] | 1050 |
|
---|
| 1051 | ///////////////////////////////////////////////////////////////
|
---|
| 1052 | #ifdef __CXX_PRAGMA_TEMPLATES__
|
---|
| 1053 | #pragma define_template LapackServer<r_4>
|
---|
| 1054 | #pragma define_template LapackServer<r_8>
|
---|
| 1055 | #pragma define_template LapackServer< complex<r_4> >
|
---|
| 1056 | #pragma define_template LapackServer< complex<r_8> >
|
---|
| 1057 | #endif
|
---|
| 1058 |
|
---|
| 1059 | #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
|
---|
| 1060 | template class LapackServer<r_4>;
|
---|
| 1061 | template class LapackServer<r_8>;
|
---|
| 1062 | template class LapackServer< complex<r_4> >;
|
---|
| 1063 | template class LapackServer< complex<r_8> >;
|
---|
| 1064 | #endif
|
---|
| 1065 |
|
---|
| 1066 | #if defined(OS_LINUX)
|
---|
| 1067 | // Pour le link avec f2c sous Linux
|
---|
| 1068 | extern "C" {
|
---|
| 1069 | void MAIN__();
|
---|
| 1070 | }
|
---|
| 1071 |
|
---|
| 1072 | void MAIN__()
|
---|
| 1073 | {
|
---|
| 1074 | cerr << "MAIN__() function for linking with libf2c.a " << endl;
|
---|
| 1075 | cerr << " This function should never be called !!! " << endl;
|
---|
| 1076 | throw PError("MAIN__() should not be called - see intflapack.cc");
|
---|
| 1077 | }
|
---|
| 1078 | #endif
|
---|