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