source: Sophya/trunk/SophyaExt/LinAlg/intflapack.cc@ 2615

Last change on this file since 2615 was 2615, checked in by cmv, 21 years ago

using namespace sophya enleve de machdefs.h, nouveau sopnamsp.h cmv 10/09/2004

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