1 | #include <iostream>
|
---|
2 | #include "intflapack.h"
|
---|
3 | #include "tvector.h"
|
---|
4 | #include "tmatrix.h"
|
---|
5 | #include <typeinfo>
|
---|
6 |
|
---|
7 | /*!
|
---|
8 | \defgroup LinAlg LinAlg module
|
---|
9 | This module contains classes and functions for complex linear
|
---|
10 | algebra on arrays. This module is intended mainly to have
|
---|
11 | classes implementing C++ interfaces between Sophya objects
|
---|
12 | and external linear algebra libraries, such as LAPACK.
|
---|
13 | */
|
---|
14 |
|
---|
15 | /*!
|
---|
16 | \class SOPHYA::LapackServer
|
---|
17 | \ingroup LinAlg
|
---|
18 | This class implements an interface to LAPACK library driver routines.
|
---|
19 | The LAPACK (Linear Algebra PACKage) is a collection high performance
|
---|
20 | routines to solve common problems in numerical linear algebra.
|
---|
21 | its is available from http://www.netlib.org.
|
---|
22 |
|
---|
23 | The present version of our LapackServer (Feb 2001) provides only
|
---|
24 | interfaces for the linear system solver and singular value
|
---|
25 | decomposition (SVD). Only arrays with BaseArray::FortranMemoryMapping
|
---|
26 | can be handled by LapackServer. LapackServer can be instanciated
|
---|
27 | for simple and double precision real or complex array types.
|
---|
28 |
|
---|
29 | The example below shows solving a linear system A*X = B
|
---|
30 |
|
---|
31 | \code
|
---|
32 | #include "intflapack.h"
|
---|
33 | // ...
|
---|
34 | // Use FortranMemoryMapping as default
|
---|
35 | BaseArray::SetDefaultMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
36 | // Create an fill the arrays A and B
|
---|
37 | int n = 20;
|
---|
38 | Matrix A(n, n);
|
---|
39 | A = RandomSequence();
|
---|
40 | Vector X(n),B(n);
|
---|
41 | X = RandomSequence();
|
---|
42 | B = A*X;
|
---|
43 | // Solve the linear system A*X = B
|
---|
44 | LapackServer<r_8> lps;
|
---|
45 | lps.LinSolve(A,B);
|
---|
46 | // We get the result in B, which should be equal to X ...
|
---|
47 | // Compute the difference B-X ;
|
---|
48 | Vector diff = B-X;
|
---|
49 | \endcode
|
---|
50 |
|
---|
51 | */
|
---|
52 |
|
---|
53 | extern "C" {
|
---|
54 | // Le calculateur de workingspace
|
---|
55 | int_4 ilaenv_(int_4 *ispec,char *name,char *opts,int_4 *n1,int_4 *n2,int_4 *n3,int_4 *n4,
|
---|
56 | int_4 nc1,int_4 nc2);
|
---|
57 |
|
---|
58 | // Drivers pour resolution de systemes lineaires
|
---|
59 | void sgesv_(int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
60 | int_4* ipiv, r_4* b, int_4* ldb, int_4* info);
|
---|
61 | void dgesv_(int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
62 | int_4* ipiv, r_8* b, int_4* ldb, int_4* info);
|
---|
63 | void cgesv_(int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
64 | int_4* ipiv, complex<r_4>* b, int_4* ldb, int_4* info);
|
---|
65 | void zgesv_(int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
66 | int_4* ipiv, complex<r_8>* b, int_4* ldb, int_4* info);
|
---|
67 |
|
---|
68 | // Drivers pour resolution de systemes lineaires symetriques
|
---|
69 | void ssysv_(char* uplo, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
70 | int_4* ipiv, r_4* b, int_4* ldb,
|
---|
71 | r_4* work, int_4* lwork, int_4* info);
|
---|
72 | void dsysv_(char* uplo, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
73 | int_4* ipiv, r_8* b, int_4* ldb,
|
---|
74 | r_8* work, int_4* lwork, int_4* info);
|
---|
75 | void csysv_(char* uplo, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
76 | int_4* ipiv, complex<r_4>* b, int_4* ldb,
|
---|
77 | complex<r_4>* work, int_4* lwork, int_4* info);
|
---|
78 | void zsysv_(char* uplo, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
79 | int_4* ipiv, complex<r_8>* b, int_4* ldb,
|
---|
80 | complex<r_8>* work, int_4* lwork, int_4* info);
|
---|
81 |
|
---|
82 | // Driver pour resolution de systemes au sens de Xi2
|
---|
83 | void sgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, r_4* a, int_4* lda,
|
---|
84 | r_4* b, int_4* ldb, r_4* work, int_4* lwork, int_4* info);
|
---|
85 | void dgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, r_8* a, int_4* lda,
|
---|
86 | r_8* b, int_4* ldb, r_8* work, int_4* lwork, int_4* info);
|
---|
87 | void cgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_4>* a, int_4* lda,
|
---|
88 | complex<r_4>* b, int_4* ldb, complex<r_4>* work, int_4* lwork, int_4* info);
|
---|
89 | void zgels_(char * trans, int_4* m, int_4* n, int_4* nrhs, complex<r_8>* a, int_4* lda,
|
---|
90 | complex<r_8>* b, int_4* ldb, complex<r_8>* work, int_4* lwork, int_4* info);
|
---|
91 |
|
---|
92 | // Driver pour decomposition SVD
|
---|
93 | void sgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_4* a, int_4* lda,
|
---|
94 | r_4* s, r_4* u, int_4* ldu, r_4* vt, int_4* ldvt,
|
---|
95 | r_4* work, int_4* lwork, int_4* info);
|
---|
96 | void dgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_8* a, int_4* lda,
|
---|
97 | r_8* s, r_8* u, int_4* ldu, r_8* vt, int_4* ldvt,
|
---|
98 | r_8* work, int_4* lwork, int_4* info);
|
---|
99 | void cgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_4>* a, int_4* lda,
|
---|
100 | complex<r_4>* s, complex<r_4>* u, int_4* ldu, complex<r_4>* vt, int_4* ldvt,
|
---|
101 | complex<r_4>* work, int_4* lwork, int_4* info);
|
---|
102 | void zgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex<r_8>* a, int_4* lda,
|
---|
103 | complex<r_8>* s, complex<r_8>* u, int_4* ldu, complex<r_8>* vt, int_4* ldvt,
|
---|
104 | complex<r_8>* work, int_4* lwork, int_4* info);
|
---|
105 |
|
---|
106 | }
|
---|
107 |
|
---|
108 |
|
---|
109 | // -------------- Classe LapackServer<T> --------------
|
---|
110 |
|
---|
111 | template <class T>
|
---|
112 | LapackServer<T>::LapackServer()
|
---|
113 | {
|
---|
114 | SetWorkSpaceSizeFactor();
|
---|
115 | }
|
---|
116 |
|
---|
117 | template <class T>
|
---|
118 | LapackServer<T>::~LapackServer()
|
---|
119 | {
|
---|
120 | }
|
---|
121 |
|
---|
122 | // --- ATTENTION BUG PREVISIBLE (CMV) --- REZA A LIRE S.T.P.
|
---|
123 | // -> Cette connerie de Fortran/C interface
|
---|
124 | // Dans les routines fortran de lapack:
|
---|
125 | // Appel depuis le C avec:
|
---|
126 | // int_4 lwork = -1;
|
---|
127 | // SUBROUTINE SSYSV( UPLO,N,NRHS,A,LDA,IPIV,B,LDB,WORK,LWORK,INFO)
|
---|
128 | // INTEGER INFO, LDA, LDB, LWORK, N, NRHS
|
---|
129 | // LOGICAL LQUERY
|
---|
130 | // LQUERY = ( LWORK.EQ.-1 )
|
---|
131 | // ELSE IF( LWORK.LT.1 .AND. .NOT.LQUERY ) THEN
|
---|
132 | // ==> le test est bien interprete sous Linux mais pas sous OSF
|
---|
133 | // ==> Sous OSF "LWORK.EQ.-1" est FALSE quand on passe lwork=-1 par argument
|
---|
134 | // ==> POUR REZA: confusion entier 4 / 8 bits ??? (bizarre on l'aurait vu avant?)
|
---|
135 | template <class T>
|
---|
136 | 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)
|
---|
137 | {
|
---|
138 | int_4 nc1 = strlen(name);
|
---|
139 | int_4 nc2 = strlen(opts);
|
---|
140 | int_4 rc=0;
|
---|
141 | rc = ilaenv_(&ispec,name,opts,&n1,&n2,&n3,&n4,nc1,nc2);
|
---|
142 | //cout<<"ilaenv_en_C("<<ispec<<","<<name<<"("<<nc1<<"),"<<opts<<"("<<nc2<<"),"
|
---|
143 | // <<n1<<","<<n2<<","<<n3<<","<<n4<<") = "<<rc<<endl;
|
---|
144 | return rc;
|
---|
145 | }
|
---|
146 |
|
---|
147 | //! Interface to Lapack linear system solver driver s/d/c/zgesvd().
|
---|
148 | /*! Solve the linear system a * x = b. Input arrays
|
---|
149 | should have FortranMemory mapping (column packed).
|
---|
150 | \param a : input matrix, overwritten on output
|
---|
151 | \param b : input-output, input vector b, contains x on exit
|
---|
152 | \return : return code from lapack driver _gesv()
|
---|
153 | */
|
---|
154 | template <class T>
|
---|
155 | int LapackServer<T>::LinSolve(TArray<T>& a, TArray<T> & b)
|
---|
156 | {
|
---|
157 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
158 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b NbDimensions() != 2"));
|
---|
159 |
|
---|
160 | int_4 rowa = a.RowsKA();
|
---|
161 | int_4 cola = a.ColsKA();
|
---|
162 | int_4 rowb = b.RowsKA();
|
---|
163 | int_4 colb = b.ColsKA();
|
---|
164 | if ( a.Size(rowa) != a.Size(cola))
|
---|
165 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Not a square Array"));
|
---|
166 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
167 | throw(SzMismatchError("LapackServer::LinSolve(a,b) RowSize(a <> b) "));
|
---|
168 |
|
---|
169 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
170 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b Not Column Packed"));
|
---|
171 |
|
---|
172 | int_4 n = a.Size(rowa);
|
---|
173 | int_4 nrhs = b.Size(colb);
|
---|
174 | int_4 lda = a.Step(cola);
|
---|
175 | int_4 ldb = b.Step(colb);
|
---|
176 | int_4 info;
|
---|
177 | int_4* ipiv = new int_4[n];
|
---|
178 |
|
---|
179 | if (typeid(T) == typeid(r_4) )
|
---|
180 | sgesv_(&n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb, &info);
|
---|
181 | else if (typeid(T) == typeid(r_8) )
|
---|
182 | dgesv_(&n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb, &info);
|
---|
183 | else if (typeid(T) == typeid(complex<r_4>) )
|
---|
184 | cgesv_(&n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
---|
185 | (complex<r_4> *)b.Data(), &ldb, &info);
|
---|
186 | else if (typeid(T) == typeid(complex<r_8>) )
|
---|
187 | zgesv_(&n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
---|
188 | (complex<r_8> *)b.Data(), &ldb, &info);
|
---|
189 | else {
|
---|
190 | delete[] ipiv;
|
---|
191 | string tn = typeid(T).name();
|
---|
192 | cerr << " LapackServer::LinSolve(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
193 | throw TypeMismatchExc("LapackServer::LinSolve(a,b) - Unsupported DataType (T)");
|
---|
194 | }
|
---|
195 | delete[] ipiv;
|
---|
196 | return(info);
|
---|
197 | }
|
---|
198 |
|
---|
199 | //! Interface to Lapack linear system solver driver s/d/c/zsysvd().
|
---|
200 | /*! Solve the linear system a * x = b with a symetric. Input arrays
|
---|
201 | should have FortranMemory mapping (column packed).
|
---|
202 | \param a : input matrix symetric , overwritten on output
|
---|
203 | \param b : input-output, input vector b, contains x on exit
|
---|
204 | \return : return code from lapack driver _gesv()
|
---|
205 | */
|
---|
206 | template <class T>
|
---|
207 | int LapackServer<T>::LinSolveSym(TArray<T>& a, TArray<T> & b)
|
---|
208 | // --- REMARQUES DE CMV ---
|
---|
209 | // 1./ contrairement a ce qui est dit dans la doc, il s'agit
|
---|
210 | // de matrices SYMETRIQUES complexes et non HERMITIENNES !!!
|
---|
211 | // 2./ pourquoi les routines de LinSolve pour des matrices symetriques
|
---|
212 | // sont plus de deux fois plus lentes que les LinSolve generales ???
|
---|
213 | {
|
---|
214 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
215 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Or b NbDimensions() != 2"));
|
---|
216 | int_4 rowa = a.RowsKA();
|
---|
217 | int_4 cola = a.ColsKA();
|
---|
218 | int_4 rowb = b.RowsKA();
|
---|
219 | int_4 colb = b.ColsKA();
|
---|
220 | if ( a.Size(rowa) != a.Size(cola))
|
---|
221 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Not a square Array"));
|
---|
222 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
223 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) RowSize(a <> b) "));
|
---|
224 |
|
---|
225 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
226 | throw(SzMismatchError("LapackServer::LinSolveSym(a,b) a Or b Not Column Packed"));
|
---|
227 |
|
---|
228 | int_4 n = a.Size(rowa);
|
---|
229 | int_4 nrhs = b.Size(colb);
|
---|
230 | int_4 lda = a.Step(cola);
|
---|
231 | int_4 ldb = b.Step(colb);
|
---|
232 | int_4 info = 0;
|
---|
233 | int_4* ipiv = new int_4[n];
|
---|
234 | int_4 lwork = -1;
|
---|
235 | T * work = NULL;
|
---|
236 |
|
---|
237 | char uplo = 'U'; // char uplo = 'L';
|
---|
238 | char struplo[5]; struplo[0] = uplo; struplo[1] = '\0';
|
---|
239 |
|
---|
240 | if (typeid(T) == typeid(r_4) ) {
|
---|
241 | lwork = ilaenv_en_C(1,"SSYTRF",struplo,n,-1,-1,-1) * n;
|
---|
242 | work = new T[lwork+5];
|
---|
243 | ssysv_(&uplo, &n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb,
|
---|
244 | (r_4 *)work, &lwork, &info);
|
---|
245 | } else if (typeid(T) == typeid(r_8) ) {
|
---|
246 | lwork = ilaenv_en_C(1,"DSYTRF",struplo,n,-1,-1,-1) * n;
|
---|
247 | work = new T[lwork+5];
|
---|
248 | dsysv_(&uplo, &n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb,
|
---|
249 | (r_8 *)work, &lwork, &info);
|
---|
250 | } else if (typeid(T) == typeid(complex<r_4>) ) {
|
---|
251 | lwork = ilaenv_en_C(1,"CSYTRF",struplo,n,-1,-1,-1) * n;
|
---|
252 | work = new T[lwork+5];
|
---|
253 | csysv_(&uplo, &n, &nrhs, (complex<r_4> *)a.Data(), &lda, ipiv,
|
---|
254 | (complex<r_4> *)b.Data(), &ldb,
|
---|
255 | (complex<r_4> *)work, &lwork, &info);
|
---|
256 | } else if (typeid(T) == typeid(complex<r_8>) ) {
|
---|
257 | lwork = ilaenv_en_C(1,"ZSYTRF",struplo,n,-1,-1,-1) * n;
|
---|
258 | work = new T[lwork+5];
|
---|
259 | zsysv_(&uplo, &n, &nrhs, (complex<r_8> *)a.Data(), &lda, ipiv,
|
---|
260 | (complex<r_8> *)b.Data(), &ldb,
|
---|
261 | (complex<r_8> *)work, &lwork, &info);
|
---|
262 | } else {
|
---|
263 | delete[] work;
|
---|
264 | delete[] ipiv;
|
---|
265 | string tn = typeid(T).name();
|
---|
266 | cerr << " LapackServer::LinSolveSym(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
267 | throw TypeMismatchExc("LapackServer::LinSolveSym(a,b) - Unsupported DataType (T)");
|
---|
268 | }
|
---|
269 | delete[] work;
|
---|
270 | delete[] ipiv;
|
---|
271 | return(info);
|
---|
272 | }
|
---|
273 |
|
---|
274 | //! Interface to Lapack least squares solver driver s/d/c/zgels().
|
---|
275 | /*! Solves the linear least squares problem defined by an m-by-n matrix
|
---|
276 | \b a and an m element vector \b b .
|
---|
277 | A solution \b x to the overdetermined system of linear equations
|
---|
278 | b = a * x is computed, minimizing the norm of b-a*x.
|
---|
279 | Underdetermined systems (m<n) are not yet handled.
|
---|
280 | Inout arrays should have FortranMemory mapping (column packed).
|
---|
281 | \param a : input matrix, overwritten on output
|
---|
282 | \param b : input-output, input vector b, contains x on exit.
|
---|
283 | \return : return code from lapack driver _gels()
|
---|
284 | \warning : b is not resized.
|
---|
285 | */
|
---|
286 | /*
|
---|
287 | $CHECK$ - A faire - cas m<n
|
---|
288 | If the linear system is underdetermined, the minimum norm
|
---|
289 | solution is computed.
|
---|
290 | */
|
---|
291 |
|
---|
292 | template <class T>
|
---|
293 | int LapackServer<T>::LeastSquareSolve(TArray<T>& a, TArray<T> & b)
|
---|
294 | {
|
---|
295 | if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) )
|
---|
296 | throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b NbDimensions() != 2"));
|
---|
297 |
|
---|
298 | int_4 rowa = a.RowsKA();
|
---|
299 | int_4 cola = a.ColsKA();
|
---|
300 | int_4 rowb = b.RowsKA();
|
---|
301 | int_4 colb = b.ColsKA();
|
---|
302 |
|
---|
303 |
|
---|
304 | if ( a.Size(rowa) != b.Size(rowb))
|
---|
305 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) RowSize(a <> b) "));
|
---|
306 |
|
---|
307 | if (!a.IsPacked(rowa) || !b.IsPacked(rowb))
|
---|
308 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) a Or b Not Column Packed"));
|
---|
309 |
|
---|
310 | if ( a.Size(rowa) < a.Size(cola)) { // $CHECK$ - m<n a changer
|
---|
311 | cout << " LapackServer<T>::LeastSquareSolve() - m<n - Not yet implemented for "
|
---|
312 | << " underdetermined systems ! " << endl;
|
---|
313 | throw(SzMismatchError("LapackServer::LeastSquareSolve(a,b) NRows<NCols - "));
|
---|
314 | }
|
---|
315 | int_4 m = a.Size(rowa);
|
---|
316 | int_4 n = a.Size(cola);
|
---|
317 | int_4 nrhs = b.Size(colb);
|
---|
318 |
|
---|
319 | int_4 lda = a.Step(cola);
|
---|
320 | int_4 ldb = b.Step(colb);
|
---|
321 | int_4 info;
|
---|
322 |
|
---|
323 | int_4 minmn = (m < n) ? m : n;
|
---|
324 | int_4 maxmn = (m > n) ? m : n;
|
---|
325 | int_4 maxmnrhs = (nrhs > maxmn) ? nrhs : maxmn;
|
---|
326 | if (maxmnrhs < 1) maxmnrhs = 1;
|
---|
327 |
|
---|
328 | int_4 lwork = minmn+maxmnrhs*5;
|
---|
329 | T * work = new T[lwork];
|
---|
330 |
|
---|
331 | char trans = 'N';
|
---|
332 |
|
---|
333 | if (typeid(T) == typeid(r_4) )
|
---|
334 | sgels_(&trans, &m, &n, &nrhs, (r_4 *)a.Data(), &lda,
|
---|
335 | (r_4 *)b.Data(), &ldb, (r_4 *)work, &lwork, &info);
|
---|
336 | else if (typeid(T) == typeid(r_8) )
|
---|
337 | dgels_(&trans, &m, &n, &nrhs, (r_8 *)a.Data(), &lda,
|
---|
338 | (r_8 *)b.Data(), &ldb, (r_8 *)work, &lwork, &info);
|
---|
339 | else if (typeid(T) == typeid(complex<r_4>) )
|
---|
340 | cgels_(&trans, &m, &n, &nrhs, (complex<r_4> *)a.Data(), &lda,
|
---|
341 | (complex<r_4> *)b.Data(), &ldb, (complex<r_4> *)work, &lwork, &info);
|
---|
342 | else if (typeid(T) == typeid(complex<r_8>) )
|
---|
343 | zgels_(&trans, &m, &n, &nrhs, (complex<r_8> *)a.Data(), &lda,
|
---|
344 | (complex<r_8> *)b.Data(), &ldb, (complex<r_8> *)work, &lwork, &info);
|
---|
345 | else {
|
---|
346 | delete[] work;
|
---|
347 | string tn = typeid(T).name();
|
---|
348 | cerr << " LapackServer::LeastSquareSolve(a,b) - Unsupported DataType T = " << tn << endl;
|
---|
349 | throw TypeMismatchExc("LapackServer::LeastSquareSolve(a,b) - Unsupported DataType (T)");
|
---|
350 | }
|
---|
351 | delete[] work;
|
---|
352 | return(info);
|
---|
353 | }
|
---|
354 |
|
---|
355 |
|
---|
356 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
357 | /*! Computes the vector of singular values of \b a. Input arrays
|
---|
358 | should have FortranMemoryMapping (column packed).
|
---|
359 | \param a : input m-by-n matrix
|
---|
360 | \param s : Vector of min(m,n) singular values (descending order)
|
---|
361 | \return : return code from lapack driver _gesvd()
|
---|
362 | */
|
---|
363 |
|
---|
364 | template <class T>
|
---|
365 | int LapackServer<T>::SVD(TArray<T>& a, TArray<T> & s)
|
---|
366 | {
|
---|
367 | return (SVDDriver(a, s, NULL, NULL) );
|
---|
368 | }
|
---|
369 |
|
---|
370 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
371 | /*! Computes the vector of singular values of \b a, as well as
|
---|
372 | right and left singular vectors of \b a.
|
---|
373 | \f[
|
---|
374 | A = U \Sigma V^T , ( A = U \Sigma V^H \ complex)
|
---|
375 | \f]
|
---|
376 | \f[
|
---|
377 | A v_i = \sigma_i u_i \ and A^T u_i = \sigma_i v_i \ (A^H \ complex)
|
---|
378 | \f]
|
---|
379 | U and V are orthogonal (unitary) matrices.
|
---|
380 | \param a : input m-by-n matrix (in FotranMemoryMapping)
|
---|
381 | \param s : Vector of min(m,n) singular values (descending order)
|
---|
382 | \param u : Matrix of left singular vectors
|
---|
383 | \param vt : Transpose of right singular vectors.
|
---|
384 | \return : return code from lapack driver _gesvd()
|
---|
385 | */
|
---|
386 | template <class T>
|
---|
387 | int LapackServer<T>::SVD(TArray<T>& a, TArray<T> & s, TArray<T> & u, TArray<T> & vt)
|
---|
388 | {
|
---|
389 | return (SVDDriver(a, s, &u, &vt) );
|
---|
390 | }
|
---|
391 |
|
---|
392 |
|
---|
393 | //! Interface to Lapack SVD driver s/d/c/zgesv().
|
---|
394 | template <class T>
|
---|
395 | int LapackServer<T>::SVDDriver(TArray<T>& a, TArray<T> & s, TArray<T>* up, TArray<T>* vtp)
|
---|
396 | {
|
---|
397 | if ( ( a.NbDimensions() != 2 ) )
|
---|
398 | throw(SzMismatchError("LapackServer::SVD(a, ...) a.NbDimensions() != 2"));
|
---|
399 |
|
---|
400 | int_4 rowa = a.RowsKA();
|
---|
401 | int_4 cola = a.ColsKA();
|
---|
402 |
|
---|
403 | if ( !a.IsPacked(rowa) )
|
---|
404 | throw(SzMismatchError("LapackServer::SVD(a, ...) a Not Column Packed "));
|
---|
405 |
|
---|
406 | int_4 m = a.Size(rowa);
|
---|
407 | int_4 n = a.Size(cola);
|
---|
408 | int_4 maxmn = (m > n) ? m : n;
|
---|
409 | int_4 minmn = (m < n) ? m : n;
|
---|
410 |
|
---|
411 | char jobu, jobvt;
|
---|
412 | jobu = 'N';
|
---|
413 | jobvt = 'N';
|
---|
414 |
|
---|
415 | sa_size_t sz[2];
|
---|
416 | if ( up != NULL) {
|
---|
417 | if ( dynamic_cast< TVector<T> * > (vtp) )
|
---|
418 | throw( TypeMismatchExc("LapackServer::SVD() Wrong type (=TVector<T>) for u !") );
|
---|
419 | up->SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
420 | sz[0] = sz[1] = m;
|
---|
421 | up->ReSize(2, sz );
|
---|
422 | jobu = 'A';
|
---|
423 | }
|
---|
424 | else {
|
---|
425 | up = new TMatrix<T>(1,1);
|
---|
426 | jobu = 'N';
|
---|
427 | }
|
---|
428 | if ( vtp != NULL) {
|
---|
429 | if ( dynamic_cast< TVector<T> * > (vtp) )
|
---|
430 | throw( TypeMismatchExc("LapackServer::SVD() Wrong type (=TVector<T>) for vt !") );
|
---|
431 | vtp->SetMemoryMapping(BaseArray::FortranMemoryMapping);
|
---|
432 | sz[0] = sz[1] = n;
|
---|
433 | vtp->ReSize(2, sz );
|
---|
434 | jobvt = 'A';
|
---|
435 | }
|
---|
436 | else {
|
---|
437 | vtp = new TMatrix<T>(1,1);
|
---|
438 | jobvt = 'N';
|
---|
439 | }
|
---|
440 |
|
---|
441 | TVector<T> *vs = dynamic_cast< TVector<T> * > (&s);
|
---|
442 | if (vs) vs->ReSize(minmn);
|
---|
443 | else {
|
---|
444 | TMatrix<T> *ms = dynamic_cast< TMatrix<T> * > (&s);
|
---|
445 | if (ms) ms->ReSize(minmn,1);
|
---|
446 | else {
|
---|
447 | sz[0] = minmn; sz[1] = 1;
|
---|
448 | s.ReSize(1, sz);
|
---|
449 | }
|
---|
450 | }
|
---|
451 |
|
---|
452 | int_4 lda = a.Step(a.ColsKA());
|
---|
453 | int_4 ldu = up->Step(up->ColsKA());
|
---|
454 | int_4 ldvt = vtp->Step(vtp->ColsKA());
|
---|
455 |
|
---|
456 | int_4 lwork = maxmn*5*wspace_size_factor;
|
---|
457 | T * work = new T[lwork];
|
---|
458 | int_4 info;
|
---|
459 |
|
---|
460 | if (typeid(T) == typeid(r_4) )
|
---|
461 | sgesvd_(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda,
|
---|
462 | (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt,
|
---|
463 | (r_4 *)work, &lwork, &info);
|
---|
464 | else if (typeid(T) == typeid(r_8) )
|
---|
465 | dgesvd_(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda,
|
---|
466 | (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt,
|
---|
467 | (r_8 *)work, &lwork, &info);
|
---|
468 | else if (typeid(T) == typeid(complex<r_4>) )
|
---|
469 | cgesvd_(&jobu, &jobvt, &m, &n, (complex<r_4> *)a.Data(), &lda,
|
---|
470 | (complex<r_4> *)s.Data(), (complex<r_4> *) up->Data(), &ldu,
|
---|
471 | (complex<r_4> *)vtp->Data(), &ldvt,
|
---|
472 | (complex<r_4> *)work, &lwork, &info);
|
---|
473 | else if (typeid(T) == typeid(complex<r_8>) )
|
---|
474 | zgesvd_(&jobu, &jobvt, &m, &n, (complex<r_8> *)a.Data(), &lda,
|
---|
475 | (complex<r_8> *)s.Data(), (complex<r_8> *) up->Data(), &ldu,
|
---|
476 | (complex<r_8> *)vtp->Data(), &ldvt,
|
---|
477 | (complex<r_8> *)work, &lwork, &info);
|
---|
478 | else {
|
---|
479 | if (jobu == 'N') delete up;
|
---|
480 | if (jobvt == 'N') delete vtp;
|
---|
481 | string tn = typeid(T).name();
|
---|
482 | cerr << " LapackServer::SVDDriver(...) - Unsupported DataType T = " << tn << endl;
|
---|
483 | throw TypeMismatchExc("LapackServer::LinSolve(a,b) - Unsupported DataType (T)");
|
---|
484 | }
|
---|
485 |
|
---|
486 | if (jobu == 'N') delete up;
|
---|
487 | if (jobvt == 'N') delete vtp;
|
---|
488 | return(info);
|
---|
489 | }
|
---|
490 |
|
---|
491 | void rztest_lapack(TArray<r_4>& aa, TArray<r_4>& bb)
|
---|
492 | {
|
---|
493 | if ( aa.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a Matrix"));
|
---|
494 | if ( aa.SizeX() != aa.SizeY()) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a square Matrix"));
|
---|
495 | if ( bb.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A Not a Matrix"));
|
---|
496 | if ( bb.SizeX() != aa.SizeX() ) throw(SzMismatchError("rztest_lapack(TMatrix<r_4> A <> B "));
|
---|
497 | if ( !bb.IsPacked() || !bb.IsPacked() )
|
---|
498 | throw(SzMismatchError("rztest_lapack(TMatrix<r_4> Not packed A or B "));
|
---|
499 |
|
---|
500 | int_4 n = aa.SizeX();
|
---|
501 | int_4 nrhs = bb.SizeY();
|
---|
502 | int_4 lda = n;
|
---|
503 | int_4 ldb = bb.SizeX();
|
---|
504 | int_4 info;
|
---|
505 | int_4* ipiv = new int_4[n];
|
---|
506 | sgesv_(&n, &nrhs, aa.Data(), &lda, ipiv, bb.Data(), &ldb, &info);
|
---|
507 | delete[] ipiv;
|
---|
508 | cout << "rztest_lapack/Info= " << info << endl;
|
---|
509 | cout << aa << "\n" << bb << endl;
|
---|
510 | return;
|
---|
511 | }
|
---|
512 |
|
---|
513 | ///////////////////////////////////////////////////////////////
|
---|
514 | #ifdef __CXX_PRAGMA_TEMPLATES__
|
---|
515 | #pragma define_template LapackServer<r_4>
|
---|
516 | #pragma define_template LapackServer<r_8>
|
---|
517 | #pragma define_template LapackServer< complex<r_4> >
|
---|
518 | #pragma define_template LapackServer< complex<r_8> >
|
---|
519 | #endif
|
---|
520 |
|
---|
521 | #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
|
---|
522 | template class LapackServer<r_4>;
|
---|
523 | template class LapackServer<r_8>;
|
---|
524 | template class LapackServer< complex<r_4> >;
|
---|
525 | template class LapackServer< complex<r_8> >;
|
---|
526 | #endif
|
---|
527 |
|
---|
528 | #if defined(OS_LINUX)
|
---|
529 | // Pour le link avec f2c sous Linux
|
---|
530 | extern "C" {
|
---|
531 | void MAIN__();
|
---|
532 | }
|
---|
533 |
|
---|
534 | void MAIN__()
|
---|
535 | {
|
---|
536 | cerr << "MAIN__() function for linking with libf2c.a " << endl;
|
---|
537 | cerr << " This function should never be called !!! " << endl;
|
---|
538 | throw PError("MAIN__() should not be called - see intflapack.cc");
|
---|
539 | }
|
---|
540 | #endif
|
---|