| 1 | #include <iostream>
|
|---|
| 2 | #include <math.h>
|
|---|
| 3 | #include "sopnamsp.h"
|
|---|
| 4 | #include "intflapack.h"
|
|---|
| 5 | #include "tvector.h"
|
|---|
| 6 | #include "tmatrix.h"
|
|---|
| 7 | #include <typeinfo>
|
|---|
| 8 |
|
|---|
| 9 | #define GARDMEM 5
|
|---|
| 10 |
|
|---|
| 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 |
|
|---|
| 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 |
|
|---|
| 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
|
|---|
| 55 | can be handled by LapackServer. LapackServer can be instanciated
|
|---|
| 56 | for simple and double precision real or complex array types.
|
|---|
| 57 | \warning The input array is overwritten in most cases.
|
|---|
| 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 |
|
|---|
| 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
|
|---|
| 172 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 173 | extern "C" {
|
|---|
| 174 | // Le calculateur de workingspace
|
|---|
| 175 | int_4 ilaenv(int_4 *ispec,char *name,char *opts,int_4 *n1,int_4 *n2,int_4 *n3,int_4 *n4,
|
|---|
| 176 | int_4 nc1,int_4 nc2);
|
|---|
| 177 |
|
|---|
| 178 | // Drivers pour resolution de systemes lineaires
|
|---|
| 179 | void sgesv(int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
|---|
| 180 | int_4* ipiv, r_4* b, int_4* ldb, int_4* info);
|
|---|
| 181 | void dgesv(int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
|---|
| 182 | int_4* ipiv, r_8* b, int_4* ldb, int_4* info);
|
|---|
| 183 | void cgesv(int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
|---|
| 184 | int_4* ipiv, complex<r_4>* b, int_4* ldb, int_4* info);
|
|---|
| 185 | void zgesv(int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
|---|
| 186 | int_4* ipiv, complex<r_8>* b, int_4* ldb, int_4* info);
|
|---|
| 187 |
|
|---|
| 188 | // Drivers pour resolution de systemes lineaires symetriques
|
|---|
| 189 | void ssysv(char* uplo, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
|---|
| 190 | int_4* ipiv, r_4* b, int_4* ldb,
|
|---|
| 191 | r_4* work, int_4* lwork, int_4* info);
|
|---|
| 192 | void dsysv(char* uplo, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
|---|
| 193 | int_4* ipiv, r_8* b, int_4* ldb,
|
|---|
| 194 | r_8* work, int_4* lwork, int_4* info);
|
|---|
| 195 | void csysv(char* uplo, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
|---|
| 196 | int_4* ipiv, complex<r_4>* b, int_4* ldb,
|
|---|
| 197 | complex<r_4>* work, int_4* lwork, int_4* info);
|
|---|
| 198 | void zsysv(char* uplo, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
|---|
| 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
|
|---|
| 203 | void sgels(char * trans, int_4* m, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
|---|
| 204 | r_4* b, int_4* ldb, r_4* work, int_4* lwork, int_4* info);
|
|---|
| 205 | void dgels(char * trans, int_4* m, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
|---|
| 206 | r_8* b, int_4* ldb, r_8* work, int_4* lwork, int_4* info);
|
|---|
| 207 | void cgels(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
|---|
| 208 | complex<r_4>* b, int_4* ldb, complex<r_4>* work, int_4* lwork, int_4* info);
|
|---|
| 209 | void zgels(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
|---|
| 210 | complex<r_8>* b, int_4* ldb, complex<r_8>* work, int_4* lwork, int_4* info);
|
|---|
| 211 |
|
|---|
| 212 | // Driver pour resolution de systemes au sens de Xi2 par SVD Divide & Conquer
|
|---|
| 213 | void sgelsd(int_4* m,int_4* n,int_4* nrhs,r_4* a,int_4* lda,
|
|---|
| 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);
|
|---|
| 216 | void dgelsd(int_4* m,int_4* n,int_4* nrhs,r_8* a,int_4* lda,
|
|---|
| 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);
|
|---|
| 219 | void cgelsd(int_4* m,int_4* n,int_4* nrhs,complex<r_4>* a,int_4* lda,
|
|---|
| 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);
|
|---|
| 222 | void zgelsd(int_4* m,int_4* n,int_4* nrhs,complex<r_8>* a,int_4* lda,
|
|---|
| 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 |
|
|---|
| 226 | // Driver pour decomposition SVD
|
|---|
| 227 | void sgesvd(char* jobu, char* jobvt, int_4* m, int_4* n, r_4* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 230 | void dgesvd(char* jobu, char* jobvt, int_4* m, int_4* n, r_8* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 233 | void cgesvd(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_4>* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 236 | void zgesvd(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_8>* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 239 |
|
|---|
| 240 | // Driver pour decomposition SVD Divide and Conquer
|
|---|
| 241 | void sgesdd(char* jobz, int_4* m, int_4* n, r_4* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 244 | void dgesdd(char* jobz, int_4* m, int_4* n, r_8* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 247 | void cgesdd(char* jobz, int_4* m, int_4* n, complex<r_4>* a, int_4* lda,
|
|---|
| 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);
|
|---|
| 250 | void zgesdd(char* jobz, int_4* m, int_4* n, complex<r_8>* a, int_4* lda,
|
|---|
| 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 |
|
|---|
| 254 | // Driver pour eigen decomposition for symetric/hermitian matrices
|
|---|
| 255 | void ssyev(char* jobz, char* uplo, int_4* n, r_4* a, int_4* lda, r_4* w,
|
|---|
| 256 | r_4* work, int_4 *lwork, int_4* info);
|
|---|
| 257 | void dsyev(char* jobz, char* uplo, int_4* n, r_8* a, int_4* lda, r_8* w,
|
|---|
| 258 | r_8* work, int_4 *lwork, int_4* info);
|
|---|
| 259 | void cheev(char* jobz, char* uplo, int_4* n, complex<r_4>* a, int_4* lda, r_4* w,
|
|---|
| 260 | complex<r_4>* work, int_4 *lwork, r_4* rwork, int_4* info);
|
|---|
| 261 | void zheev(char* jobz, char* uplo, int_4* n, complex<r_8>* a, int_4* lda, r_8* w,
|
|---|
| 262 | complex<r_8>* work, int_4 *lwork, r_8* rwork, int_4* info);
|
|---|
| 263 |
|
|---|
| 264 | // Driver pour eigen decomposition for general squared matrices
|
|---|
| 265 | void sgeev(char* jobl, char* jobvr, int_4* n, r_4* a, int_4* lda, r_4* wr, r_4* wi,
|
|---|
| 266 | r_4* vl, int_4* ldvl, r_4* vr, int_4* ldvr,
|
|---|
| 267 | r_4* work, int_4 *lwork, int_4* info);
|
|---|
| 268 | void dgeev(char* jobl, char* jobvr, int_4* n, r_8* a, int_4* lda, r_8* wr, r_8* wi,
|
|---|
| 269 | r_8* vl, int_4* ldvl, r_8* vr, int_4* ldvr,
|
|---|
| 270 | r_8* work, int_4 *lwork, int_4* info);
|
|---|
| 271 | void cgeev(char* jobl, char* jobvr, int_4* n, complex<r_4>* a, int_4* lda, complex<r_4>* w,
|
|---|
| 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);
|
|---|
| 274 | void zgeev(char* jobl, char* jobvr, int_4* n, complex<r_8>* a, int_4* lda, complex<r_8>* w,
|
|---|
| 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 |
|
|---|
| 278 | }
|
|---|
| 279 |
|
|---|
| 280 | // -------------- Classe LapackServer<T> --------------
|
|---|
| 281 |
|
|---|
| 282 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 283 | template <class T>
|
|---|
| 284 | LapackServer<T>::LapackServer()
|
|---|
| 285 | {
|
|---|
| 286 | SetWorkSpaceSizeFactor();
|
|---|
| 287 | }
|
|---|
| 288 |
|
|---|
| 289 | template <class T>
|
|---|
| 290 | LapackServer<T>::~LapackServer()
|
|---|
| 291 | {
|
|---|
| 292 | }
|
|---|
| 293 |
|
|---|
| 294 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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 | {
|
|---|
| 298 | int_4 nc1 = strlen(name), nc2 = strlen(opts), rc=0;
|
|---|
| 299 | rc = ilaenv(&ispec,name,opts,&n1,&n2,&n3,&n4,nc1,nc2);
|
|---|
| 300 | //cout<<"ilaenv_en_C("<<ispec<<","<<name<<"("<<nc1<<"),"<<opts<<"("<<nc2<<"),"
|
|---|
| 301 | // <<n1<<","<<n2<<","<<n3<<","<<n4<<") = "<<rc<<endl;
|
|---|
| 302 | return rc;
|
|---|
| 303 | }
|
|---|
| 304 |
|
|---|
| 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 |
|
|---|
| 320 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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).
|
|---|
| 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 | */
|
|---|
| 328 | template <class T>
|
|---|
| 329 | int LapackServer<T>::LinSolve(TArray<T>& a, TArray<T> & b)
|
|---|
| 330 | {
|
|---|
| 331 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
|---|
| 332 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b NbDimensions() != 2"));
|
|---|
| 333 |
|
|---|
| 334 | int_4 rowa = a.RowsKA();
|
|---|
| 335 | int_4 cola = a.ColsKA();
|
|---|
| 336 | int_4 rowb = b.RowsKA();
|
|---|
| 337 | int_4 colb = b.ColsKA();
|
|---|
| 338 | if ( a.Size(rowa) != a.Size(cola))
|
|---|
| 339 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Not a square Array"));
|
|---|
| 340 | if ( a.Size(rowa) != b.Size(rowb))
|
|---|
| 341 | throw(SzMismatchError("LapackServer::LinSolve(a,b) RowSize(a <> b) "));
|
|---|
| 342 |
|
|---|
| 343 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
|---|
| 344 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b Not Column Packed"));
|
|---|
| 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) )
|
|---|
| 354 | sgesv(&n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb, &info);
|
|---|
| 355 | else if (typeid(T) == typeid(r_8) )
|
|---|
| 356 | dgesv(&n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb, &info);
|
|---|
| 357 | else if (typeid(T) == typeid(complex<r_4>) )
|
|---|
| 358 | cgesv(&n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
|---|
| 359 | (complex<r_4> *)b.Data(), &ldb, &info);
|
|---|
| 360 | else if (typeid(T) == typeid(complex<r_8>) )
|
|---|
| 361 | zgesv(&n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
|---|
| 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;
|
|---|
| 370 | return(info);
|
|---|
| 371 | }
|
|---|
| 372 |
|
|---|
| 373 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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).
|
|---|
| 377 | \param a : input matrix symetric , overwritten on output
|
|---|
| 378 | \param b : input-output, input vector b, contains x on exit
|
|---|
| 379 | \return : return code from lapack driver _sysv()
|
|---|
| 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
|
|---|
| 387 | // sont plus de deux fois plus lentes que les LinSolve generales sur OSF
|
|---|
| 388 | // et sensiblement plus lentes sous Linux ???
|
|---|
| 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;
|
|---|
| 412 | T wkget[2];
|
|---|
| 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) ) {
|
|---|
| 418 | ssysv(&uplo, &n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb,
|
|---|
| 419 | (r_4 *)wkget, &lwork, &info);
|
|---|
| 420 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 421 | ssysv(&uplo, &n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb,
|
|---|
| 422 | (r_4 *)work, &lwork, &info);
|
|---|
| 423 | } else if (typeid(T) == typeid(r_8) ) {
|
|---|
| 424 | dsysv(&uplo, &n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb,
|
|---|
| 425 | (r_8 *)wkget, &lwork, &info);
|
|---|
| 426 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 427 | dsysv(&uplo, &n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb,
|
|---|
| 428 | (r_8 *)work, &lwork, &info);
|
|---|
| 429 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 430 | csysv(&uplo, &n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
|---|
| 431 | (complex<r_4> *)b.Data(), &ldb,
|
|---|
| 432 | (complex<r_4> *)wkget, &lwork, &info);
|
|---|
| 433 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 434 | csysv(&uplo, &n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
|---|
| 435 | (complex<r_4> *)b.Data(), &ldb,
|
|---|
| 436 | (complex<r_4> *)work, &lwork, &info);
|
|---|
| 437 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 438 | zsysv(&uplo, &n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
|---|
| 439 | (complex<r_8> *)b.Data(), &ldb,
|
|---|
| 440 | (complex<r_8> *)wkget, &lwork, &info);
|
|---|
| 441 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 442 | zsysv(&uplo, &n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
|---|
| 443 | (complex<r_8> *)b.Data(), &ldb,
|
|---|
| 444 | (complex<r_8> *)work, &lwork, &info);
|
|---|
| 445 | } else {
|
|---|
| 446 | if(work) delete[] work;
|
|---|
| 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 | }
|
|---|
| 452 | if(work) delete[] work;
|
|---|
| 453 | delete[] ipiv;
|
|---|
| 454 | return(info);
|
|---|
| 455 | }
|
|---|
| 456 |
|
|---|
| 457 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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
|
|---|
| 460 | \b a and an m element vector \b b , using QR or LQ factorization .
|
|---|
| 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 |
|
|---|
| 476 | template <class T>
|
|---|
| 477 | int LapackServer<T>::LeastSquareSolve(TArray<T>& a, TArray<T> & b)
|
|---|
| 478 | {
|
|---|
| 479 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
|---|
| 480 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) a Or b NbDimensions() != 2"));
|
|---|
| 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))
|
|---|
| 492 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) a Or b Not Column Packed"));
|
|---|
| 493 |
|
|---|
| 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 | }
|
|---|
| 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 |
|
|---|
| 512 | int_4 lwork = -1; //minmn+maxmnrhs*5;
|
|---|
| 513 | T * work = NULL;
|
|---|
| 514 | T wkget[2];
|
|---|
| 515 |
|
|---|
| 516 | char trans = 'N';
|
|---|
| 517 |
|
|---|
| 518 | if (typeid(T) == typeid(r_4) ) {
|
|---|
| 519 | sgels(&trans, &m, &n, &nrhs, (r_4 *)a.Data(), &lda,
|
|---|
| 520 | (r_4 *)b.Data(), &ldb, (r_4 *)wkget, &lwork, &info);
|
|---|
| 521 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 522 | sgels(&trans, &m, &n, &nrhs, (r_4 *)a.Data(), &lda,
|
|---|
| 523 | (r_4 *)b.Data(), &ldb, (r_4 *)work, &lwork, &info);
|
|---|
| 524 | } else if (typeid(T) == typeid(r_8) ) {
|
|---|
| 525 | dgels(&trans, &m, &n, &nrhs, (r_8 *)a.Data(), &lda,
|
|---|
| 526 | (r_8 *)b.Data(), &ldb, (r_8 *)wkget, &lwork, &info);
|
|---|
| 527 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 528 | dgels(&trans, &m, &n, &nrhs, (r_8 *)a.Data(), &lda,
|
|---|
| 529 | (r_8 *)b.Data(), &ldb, (r_8 *)work, &lwork, &info);
|
|---|
| 530 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 531 | cgels(&trans, &m, &n, &nrhs, (complex<r_4> *)a.Data(), &lda,
|
|---|
| 532 | (complex<r_4> *)b.Data(), &ldb, (complex<r_4> *)wkget, &lwork, &info);
|
|---|
| 533 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 534 | cgels(&trans, &m, &n, &nrhs, (complex<r_4> *)a.Data(), &lda,
|
|---|
| 535 | (complex<r_4> *)b.Data(), &ldb, (complex<r_4> *)work, &lwork, &info);
|
|---|
| 536 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 537 | zgels(&trans, &m, &n, &nrhs, (complex<r_8> *)a.Data(), &lda,
|
|---|
| 538 | (complex<r_8> *)b.Data(), &ldb, (complex<r_8> *)wkget, &lwork, &info);
|
|---|
| 539 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 540 | zgels(&trans, &m, &n, &nrhs, (complex<r_8> *)a.Data(), &lda,
|
|---|
| 541 | (complex<r_8> *)b.Data(), &ldb, (complex<r_8> *)work, &lwork, &info);
|
|---|
| 542 | } else {
|
|---|
| 543 | if(work) delete [] work; work=NULL;
|
|---|
| 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 | }
|
|---|
| 548 | if(work) delete [] work;
|
|---|
| 549 | return(info);
|
|---|
| 550 | }
|
|---|
| 551 |
|
|---|
| 552 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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.
|
|---|
| 593 | \param rank : output the rank of the matrix.
|
|---|
| 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"));
|
|---|
| 605 |
|
|---|
| 606 | int_4 m = a.NRows();
|
|---|
| 607 | int_4 n = a.NCols();
|
|---|
| 608 |
|
|---|
| 609 | if(b.NRows() != m)
|
|---|
| 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 |
|
|---|
| 620 | { // Use {} for automatic des-allocation of "bsave"
|
|---|
| 621 | TMatrix<T> bsave(m,nrhs); bsave.SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
|---|
| 622 | bsave = b;
|
|---|
| 623 | b.ReSize(maxmn,nrhs); b = (T) 0;
|
|---|
| 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"
|
|---|
| 626 | s.ReSize(minmn);
|
|---|
| 627 |
|
|---|
| 628 | int_4 smlsiz = 25; // Normallement ilaenv_en_C(9,...) renvoie toujours 25
|
|---|
| 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 |
|
|---|
| 637 | T * work = NULL;
|
|---|
| 638 | int_4 * iwork = NULL;
|
|---|
| 639 | int_4 lwork=-1, lrwork;
|
|---|
| 640 | T wkget[2];
|
|---|
| 641 |
|
|---|
| 642 | if(typeid(T) == typeid(r_4) ) {
|
|---|
| 643 | r_4* sloc = new r_4[minmn];
|
|---|
| 644 | r_4 srcond = rcond;
|
|---|
| 645 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
|---|
| 646 | sgelsd(&m,&n,&nrhs,(r_4*)a.Data(),&lda,
|
|---|
| 647 | (r_4*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
|---|
| 648 | (r_4*)wkget,&lwork,(int_4*)iwork,&info);
|
|---|
| 649 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 650 | sgelsd(&m,&n,&nrhs,(r_4*)a.Data(),&lda,
|
|---|
| 651 | (r_4*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
|---|
| 652 | (r_4*)work,&lwork,(int_4*)iwork,&info);
|
|---|
| 653 | for(int_4 i=0;i<minmn;i++) s(i) = sloc[i];
|
|---|
| 654 | delete [] sloc;
|
|---|
| 655 | } else if(typeid(T) == typeid(r_8) ) {
|
|---|
| 656 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
|---|
| 657 | dgelsd(&m,&n,&nrhs,(r_8*)a.Data(),&lda,
|
|---|
| 658 | (r_8*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
|---|
| 659 | (r_8*)wkget,&lwork,(int_4*)iwork,&info);
|
|---|
| 660 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 661 | dgelsd(&m,&n,&nrhs,(r_8*)a.Data(),&lda,
|
|---|
| 662 | (r_8*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
|---|
| 663 | (r_8*)work,&lwork,(int_4*)iwork,&info);
|
|---|
| 664 | } else if(typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 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;
|
|---|
| 669 | r_4* rwork = new r_4[lrwork +GARDMEM];
|
|---|
| 670 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
|---|
| 671 | r_4* sloc = new r_4[minmn];
|
|---|
| 672 | r_4 srcond = rcond;
|
|---|
| 673 | cgelsd(&m,&n,&nrhs,(complex<r_4>*)a.Data(),&lda,
|
|---|
| 674 | (complex<r_4>*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
|---|
| 675 | (complex<r_4>*)wkget,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
|---|
| 676 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 677 | cgelsd(&m,&n,&nrhs,(complex<r_4>*)a.Data(),&lda,
|
|---|
| 678 | (complex<r_4>*)b.Data(),&ldb,(r_4*)sloc,&srcond,&rank,
|
|---|
| 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];
|
|---|
| 681 | delete [] sloc; delete [] rwork;
|
|---|
| 682 | } else if(typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 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;
|
|---|
| 688 | r_8* rwork = new r_8[lrwork +GARDMEM];
|
|---|
| 689 | iwork = new int_4[3*minmn*nlvl+11*minmn +GARDMEM];
|
|---|
| 690 | zgelsd(&m,&n,&nrhs,(complex<r_8>*)a.Data(),&lda,
|
|---|
| 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];
|
|---|
| 694 | zgelsd(&m,&n,&nrhs,(complex<r_8>*)a.Data(),&lda,
|
|---|
| 695 | (complex<r_8>*)b.Data(),&ldb,(r_8*)s.Data(),&rcond,&rank,
|
|---|
| 696 | (complex<r_8>*)work,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
|---|
| 697 | delete [] rwork;
|
|---|
| 698 | } else {
|
|---|
| 699 | if(work) delete [] work; work=NULL;
|
|---|
| 700 | if(iwork) delete [] iwork; iwork=NULL;
|
|---|
| 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 |
|
|---|
| 706 | if(work) delete [] work; if(iwork) delete [] iwork;
|
|---|
| 707 | return(info);
|
|---|
| 708 | }
|
|---|
| 709 |
|
|---|
| 710 |
|
|---|
| 711 | ////////////////////////////////////////////////////////////////////////////////////
|
|---|
| 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 |
|
|---|
| 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 |
|
|---|
| 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.
|
|---|
| 736 | \param a : input m-by-n matrix (in FortranMemoryMapping)
|
|---|
| 737 | \param s : Vector of min(m,n) singular values (descending order)
|
|---|
| 738 | \param u : m-by-m Matrix of left singular vectors
|
|---|
| 739 | \param vt : Transpose of right singular vectors (n-by-n matrix).
|
|---|
| 740 | \return : return code from lapack driver _gesvd()
|
|---|
| 741 | */
|
|---|
| 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 |
|
|---|
| 748 |
|
|---|
| 749 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
|---|
| 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 ) )
|
|---|
| 754 | throw(SzMismatchError("LapackServer::SVDDriver(a, ...) a.NbDimensions() != 2"));
|
|---|
| 755 |
|
|---|
| 756 | int_4 rowa = a.RowsKA();
|
|---|
| 757 | int_4 cola = a.ColsKA();
|
|---|
| 758 |
|
|---|
| 759 | if ( !a.IsPacked(rowa) )
|
|---|
| 760 | throw(SzMismatchError("LapackServer::SVDDriver(a, ...) a Not Column Packed "));
|
|---|
| 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) )
|
|---|
| 774 | throw( TypeMismatchExc("LapackServer::SVDDriver() Wrong type (=TVector<T>) for u !") );
|
|---|
| 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) )
|
|---|
| 786 | throw( TypeMismatchExc("LapackServer::SVDDriver() Wrong type (=TVector<T>) for vt !") );
|
|---|
| 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());
|
|---|
| 811 | int_4 info;
|
|---|
| 812 |
|
|---|
| 813 | int_4 lwork = -1; // maxmn*5 *wspace_size_factor;
|
|---|
| 814 | T * work = NULL; // = new T[lwork];
|
|---|
| 815 | T wkget[2];
|
|---|
| 816 |
|
|---|
| 817 | if (typeid(T) == typeid(r_4) ) {
|
|---|
| 818 | sgesvd(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda,
|
|---|
| 819 | (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt,
|
|---|
| 820 | (r_4 *)wkget, &lwork, &info);
|
|---|
| 821 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 822 | sgesvd(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda,
|
|---|
| 823 | (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt,
|
|---|
| 824 | (r_4 *)work, &lwork, &info);
|
|---|
| 825 | } else if (typeid(T) == typeid(r_8) ) {
|
|---|
| 826 | dgesvd(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda,
|
|---|
| 827 | (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt,
|
|---|
| 828 | (r_8 *)wkget, &lwork, &info);
|
|---|
| 829 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 830 | dgesvd(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda,
|
|---|
| 831 | (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt,
|
|---|
| 832 | (r_8 *)work, &lwork, &info);
|
|---|
| 833 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 834 | r_4 * rwork = new r_4[5*minmn +GARDMEM];
|
|---|
| 835 | r_4 * sloc = new r_4[minmn];
|
|---|
| 836 | cgesvd(&jobu, &jobvt, &m, &n, (complex<r_4> *)a.Data(), &lda,
|
|---|
| 837 | (r_4 *)sloc, (complex<r_4> *) up->Data(), &ldu,
|
|---|
| 838 | (complex<r_4> *)vtp->Data(), &ldvt,
|
|---|
| 839 | (complex<r_4> *)wkget, &lwork, (r_4 *)rwork, &info);
|
|---|
| 840 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 841 | cgesvd(&jobu, &jobvt, &m, &n, (complex<r_4> *)a.Data(), &lda,
|
|---|
| 842 | (r_4 *)sloc, (complex<r_4> *) up->Data(), &ldu,
|
|---|
| 843 | (complex<r_4> *)vtp->Data(), &ldvt,
|
|---|
| 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>) ) {
|
|---|
| 848 | r_8 * rwork = new r_8[5*minmn +GARDMEM];
|
|---|
| 849 | r_8 * sloc = new r_8[minmn];
|
|---|
| 850 | zgesvd(&jobu, &jobvt, &m, &n, (complex<r_8> *)a.Data(), &lda,
|
|---|
| 851 | (r_8 *)sloc, (complex<r_8> *) up->Data(), &ldu,
|
|---|
| 852 | (complex<r_8> *)vtp->Data(), &ldvt,
|
|---|
| 853 | (complex<r_8> *)wkget, &lwork, (r_8 *)rwork, &info);
|
|---|
| 854 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 855 | zgesvd(&jobu, &jobvt, &m, &n, (complex<r_8> *)a.Data(), &lda,
|
|---|
| 856 | (r_8 *)sloc, (complex<r_8> *) up->Data(), &ldu,
|
|---|
| 857 | (complex<r_8> *)vtp->Data(), &ldvt,
|
|---|
| 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 {
|
|---|
| 862 | if(work) delete [] work; work=NULL;
|
|---|
| 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;
|
|---|
| 867 | throw TypeMismatchExc("LapackServer::SVDDriver(a,b) - Unsupported DataType (T)");
|
|---|
| 868 | }
|
|---|
| 869 |
|
|---|
| 870 | if(work) delete [] work;
|
|---|
| 871 | if (jobu == 'N') delete up;
|
|---|
| 872 | if (jobvt == 'N') delete vtp;
|
|---|
| 873 | return(info);
|
|---|
| 874 | }
|
|---|
| 875 |
|
|---|
| 876 |
|
|---|
| 877 | //! Interface to Lapack SVD driver s/d/c/zgesdd().
|
|---|
| 878 | /*! Same as SVD but with Divide and Conquer method */
|
|---|
| 879 | template <class T>
|
|---|
| 880 | int LapackServer<T>::SVD_DC(TMatrix<T>& a, TVector<r_8>& s, TMatrix<T>& u, TMatrix<T>& vt)
|
|---|
| 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 |
|
|---|
| 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 |
|
|---|
| 907 | if(typeid(T) == typeid(r_4) ) {
|
|---|
| 908 | r_4* sloc = new r_4[minmn];
|
|---|
| 909 | iwork = new int_4[8*minmn +GARDMEM];
|
|---|
| 910 | sgesdd(&jobz,&m,&n,(r_4*)a.Data(),&lda,
|
|---|
| 911 | (r_4*)sloc,(r_4*)u.Data(),&ldu,(r_4*)vt.Data(),&ldvt,
|
|---|
| 912 | (r_4*)wkget,&lwork,(int_4*)iwork,&info);
|
|---|
| 913 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 914 | sgesdd(&jobz,&m,&n,(r_4*)a.Data(),&lda,
|
|---|
| 915 | (r_4*)sloc,(r_4*)u.Data(),&ldu,(r_4*)vt.Data(),&ldvt,
|
|---|
| 916 | (r_4*)work,&lwork,(int_4*)iwork,&info);
|
|---|
| 917 | for(int_4 i=0;i<minmn;i++) s(i) = (r_8) sloc[i];
|
|---|
| 918 | delete [] sloc;
|
|---|
| 919 | } else if(typeid(T) == typeid(r_8) ) {
|
|---|
| 920 | iwork = new int_4[8*minmn +GARDMEM];
|
|---|
| 921 | dgesdd(&jobz,&m,&n,(r_8*)a.Data(),&lda,
|
|---|
| 922 | (r_8*)s.Data(),(r_8*)u.Data(),&ldu,(r_8*)vt.Data(),&ldvt,
|
|---|
| 923 | (r_8*)wkget,&lwork,(int_4*)iwork,&info);
|
|---|
| 924 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 925 | dgesdd(&jobz,&m,&n,(r_8*)a.Data(),&lda,
|
|---|
| 926 | (r_8*)s.Data(),(r_8*)u.Data(),&ldu,(r_8*)vt.Data(),&ldvt,
|
|---|
| 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];
|
|---|
| 930 | r_4* rwork = new r_4[5*minmn*minmn+5*minmn +GARDMEM];
|
|---|
| 931 | iwork = new int_4[8*minmn +GARDMEM];
|
|---|
| 932 | cgesdd(&jobz,&m,&n,(complex<r_4>*)a.Data(),&lda,
|
|---|
| 933 | (r_4*)sloc,(complex<r_4>*)u.Data(),&ldu,(complex<r_4>*)vt.Data(),&ldvt,
|
|---|
| 934 | (complex<r_4>*)wkget,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
|---|
| 935 | lwork = type2i4(&wkget[0],4); work = new T[lwork +GARDMEM];
|
|---|
| 936 | cgesdd(&jobz,&m,&n,(complex<r_4>*)a.Data(),&lda,
|
|---|
| 937 | (r_4*)sloc,(complex<r_4>*)u.Data(),&ldu,(complex<r_4>*)vt.Data(),&ldvt,
|
|---|
| 938 | (complex<r_4>*)work,&lwork,(r_4*)rwork,(int_4*)iwork,&info);
|
|---|
| 939 | for(int_4 i=0;i<minmn;i++) s(i) = (r_8) sloc[i];
|
|---|
| 940 | delete [] sloc; delete [] rwork;
|
|---|
| 941 | } else if(typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 942 | r_8* rwork = new r_8[5*minmn*minmn+5*minmn +GARDMEM];
|
|---|
| 943 | iwork = new int_4[8*minmn +GARDMEM];
|
|---|
| 944 | zgesdd(&jobz,&m,&n,(complex<r_8>*)a.Data(),&lda,
|
|---|
| 945 | (r_8*)s.Data(),(complex<r_8>*)u.Data(),&ldu,(complex<r_8>*)vt.Data(),&ldvt,
|
|---|
| 946 | (complex<r_8>*)wkget,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
|---|
| 947 | lwork = type2i4(&wkget[0],8); work = new T[lwork +GARDMEM];
|
|---|
| 948 | zgesdd(&jobz,&m,&n,(complex<r_8>*)a.Data(),&lda,
|
|---|
| 949 | (r_8*)s.Data(),(complex<r_8>*)u.Data(),&ldu,(complex<r_8>*)vt.Data(),&ldvt,
|
|---|
| 950 | (complex<r_8>*)work,&lwork,(r_8*)rwork,(int_4*)iwork,&info);
|
|---|
| 951 | delete [] rwork;
|
|---|
| 952 | } else {
|
|---|
| 953 | if(work) delete [] work; work=NULL;
|
|---|
| 954 | if(iwork) delete [] iwork; iwork=NULL;
|
|---|
| 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 |
|
|---|
| 960 | if(work) delete [] work; if(iwork) delete [] iwork;
|
|---|
| 961 | return(info);
|
|---|
| 962 | }
|
|---|
| 963 |
|
|---|
| 964 |
|
|---|
| 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)
|
|---|
| 972 | \return : return code from lapack driver
|
|---|
| 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 |
|
|---|
| 987 | char uplo='U';
|
|---|
| 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;
|
|---|
| 993 | int_4 lwork = -1;
|
|---|
| 994 | T * work = NULL;
|
|---|
| 995 | T wkget[2];
|
|---|
| 996 |
|
|---|
| 997 | b.ReSize(n); b = 0.;
|
|---|
| 998 |
|
|---|
| 999 | if (typeid(T) == typeid(r_4) ) {
|
|---|
| 1000 | r_4* w = new r_4[n];
|
|---|
| 1001 | ssyev(&jobz,&uplo,&n,(r_4 *)a.Data(),&lda,(r_4 *)w,(r_4 *)wkget,&lwork,&info);
|
|---|
| 1002 | lwork = type2i4(&wkget[0],4); /* 3*n-1;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1003 | ssyev(&jobz,&uplo,&n,(r_4 *)a.Data(),&lda,(r_4 *)w,(r_4 *)work,&lwork,&info);
|
|---|
| 1004 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
|---|
| 1005 | delete [] w;
|
|---|
| 1006 | } else if (typeid(T) == typeid(r_8) ) {
|
|---|
| 1007 | r_8* w = new r_8[n];
|
|---|
| 1008 | dsyev(&jobz,&uplo,&n,(r_8 *)a.Data(),&lda,(r_8 *)w,(r_8 *)wkget,&lwork,&info);
|
|---|
| 1009 | lwork = type2i4(&wkget[0],8); /* 3*n-1;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1010 | dsyev(&jobz,&uplo,&n,(r_8 *)a.Data(),&lda,(r_8 *)w,(r_8 *)work,&lwork,&info);
|
|---|
| 1011 | if(info==0) for(int i=0;i<n;i++) b(i) = w[i];
|
|---|
| 1012 | delete [] w;
|
|---|
| 1013 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 1014 | r_4* rwork = new r_4[3*n-2 +GARDMEM]; r_4* w = new r_4[n];
|
|---|
| 1015 | cheev(&jobz,&uplo,&n,(complex<r_4> *)a.Data(),&lda,(r_4 *)w
|
|---|
| 1016 | ,(complex<r_4> *)wkget,&lwork,(r_4 *)rwork,&info);
|
|---|
| 1017 | lwork = type2i4(&wkget[0],4); /* 2*n-1;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1018 | cheev(&jobz,&uplo,&n,(complex<r_4> *)a.Data(),&lda,(r_4 *)w
|
|---|
| 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];
|
|---|
| 1021 | delete [] rwork; delete [] w;
|
|---|
| 1022 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 1023 | r_8* rwork = new r_8[3*n-2 +GARDMEM]; r_8* w = new r_8[n];
|
|---|
| 1024 | zheev(&jobz,&uplo,&n,(complex<r_8> *)a.Data(),&lda,(r_8 *)w
|
|---|
| 1025 | ,(complex<r_8> *)wkget,&lwork,(r_8 *)rwork,&info);
|
|---|
| 1026 | lwork = type2i4(&wkget[0],8); /* 2*n-1;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1027 | zheev(&jobz,&uplo,&n,(complex<r_8> *)a.Data(),&lda,(r_8 *)w
|
|---|
| 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];
|
|---|
| 1030 | delete [] rwork; delete [] w;
|
|---|
| 1031 | } else {
|
|---|
| 1032 | if(work) delete [] work; work=NULL;
|
|---|
| 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 |
|
|---|
| 1038 | if(work) delete [] work;
|
|---|
| 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
|
|---|
| 1050 | \return : return code from lapack driver
|
|---|
| 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 |
|
|---|
| 1079 | char jobvl = 'N';
|
|---|
| 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 |
|
|---|
| 1090 | int_4 lwork = -1;
|
|---|
| 1091 | T * work = NULL;
|
|---|
| 1092 | T wkget[2];
|
|---|
| 1093 |
|
|---|
| 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;
|
|---|
| 1096 | sgeev(&jobvl,&jobvr,&n,(r_4 *)a.Data(),&lda,(r_4 *)wr,(r_4 *)wi,
|
|---|
| 1097 | (r_4 *)vl,&ldvl,(r_4 *)evec.Data(),&ldvr,
|
|---|
| 1098 | (r_4 *)wkget,&lwork,&info);
|
|---|
| 1099 | lwork = type2i4(&wkget[0],4); /* 4*n;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1100 | sgeev(&jobvl,&jobvr,&n,(r_4 *)a.Data(),&lda,(r_4 *)wr,(r_4 *)wi,
|
|---|
| 1101 | (r_4 *)vl,&ldvl,(r_4 *)evec.Data(),&ldvr,
|
|---|
| 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]);
|
|---|
| 1104 | delete [] wr; delete [] wi;
|
|---|
| 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;
|
|---|
| 1107 | dgeev(&jobvl,&jobvr,&n,(r_8 *)a.Data(),&lda,(r_8 *)wr,(r_8 *)wi,
|
|---|
| 1108 | (r_8 *)vl,&ldvl,(r_8 *)evec.Data(),&ldvr,
|
|---|
| 1109 | (r_8 *)wkget,&lwork,&info);
|
|---|
| 1110 | lwork = type2i4(&wkget[0],8); /* 4*n;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1111 | dgeev(&jobvl,&jobvr,&n,(r_8 *)a.Data(),&lda,(r_8 *)wr,(r_8 *)wi,
|
|---|
| 1112 | (r_8 *)vl,&ldvl,(r_8 *)evec.Data(),&ldvr,
|
|---|
| 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]);
|
|---|
| 1115 | delete [] wr; delete [] wi;
|
|---|
| 1116 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
|---|
| 1117 | r_4* rwork = new r_4[2*n +GARDMEM]; r_4* vl = NULL; TVector< complex<r_4> > w(n);
|
|---|
| 1118 | cgeev(&jobvl,&jobvr,&n,(complex<r_4> *)a.Data(),&lda,(complex<r_4> *)w.Data(),
|
|---|
| 1119 | (complex<r_4> *)vl,&ldvl,(complex<r_4> *)evec.Data(),&ldvr,
|
|---|
| 1120 | (complex<r_4> *)wkget,&lwork,(r_4 *)rwork,&info);
|
|---|
| 1121 | lwork = type2i4(&wkget[0],4); /* 2*n;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1122 | cgeev(&jobvl,&jobvr,&n,(complex<r_4> *)a.Data(),&lda,(complex<r_4> *)w.Data(),
|
|---|
| 1123 | (complex<r_4> *)vl,&ldvl,(complex<r_4> *)evec.Data(),&ldvr,
|
|---|
| 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);
|
|---|
| 1126 | delete [] rwork;
|
|---|
| 1127 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
|---|
| 1128 | r_8* rwork = new r_8[2*n +GARDMEM]; r_8* vl = NULL;
|
|---|
| 1129 | zgeev(&jobvl,&jobvr,&n,(complex<r_8> *)a.Data(),&lda,(complex<r_8> *)eval.Data(),
|
|---|
| 1130 | (complex<r_8> *)vl,&ldvl,(complex<r_8> *)evec.Data(),&ldvr,
|
|---|
| 1131 | (complex<r_8> *)wkget,&lwork,(r_8 *)rwork,&info);
|
|---|
| 1132 | lwork = type2i4(&wkget[0],8); /* 2*n;*/ work = new T[lwork +GARDMEM];
|
|---|
| 1133 | zgeev(&jobvl,&jobvr,&n,(complex<r_8> *)a.Data(),&lda,(complex<r_8> *)eval.Data(),
|
|---|
| 1134 | (complex<r_8> *)vl,&ldvl,(complex<r_8> *)evec.Data(),&ldvr,
|
|---|
| 1135 | (complex<r_8> *)work,&lwork,(r_8 *)rwork,&info);
|
|---|
| 1136 | delete [] rwork;
|
|---|
| 1137 | } else {
|
|---|
| 1138 | if(work) delete [] work; work=NULL;
|
|---|
| 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 |
|
|---|
| 1144 | if(work) delete [] work;
|
|---|
| 1145 | return(info);
|
|---|
| 1146 | }
|
|---|
| 1147 |
|
|---|
| 1148 |
|
|---|
| 1149 |
|
|---|
| 1150 |
|
|---|
| 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)
|
|---|
| 1160 | namespace SOPHYA {
|
|---|
| 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> >;
|
|---|
| 1165 | }
|
|---|
| 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
|
|---|