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