source: Sophya/trunk/SophyaProg/Tests/tsttminv.cc@ 2557

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

cosmetique... cmv 21/07/04

File size: 17.8 KB
RevLine 
[2555]1// Test de l'inversion de matrices et valeurs propres (avec Lapack) (cmv 21/07/04)
2// cmvtminv -a 1234 -l 0 -s 1 -b 25,10000 -n 50 -S
[995]3
[2555]4///////////////////////////////////////////////////
5///////////////////////////////////////////////////
6// PARTIE POUVANT ETRE CHANGEE PAR L'UTILISATEUR //
7///////////////////////////////////////////////////
8///////////////////////////////////////////////////
9
10// --- Choix de travailler avec des matrices complexes ?
[995]11//#define COMPLEX
12
[2555]13//////////////////////////////////////////////////
14// --- Choix de travailler en simple precision ?
15//#define PRECIS32
16
17//////////////////////////////////////////////////
18// --- Choix GausPiv + Lapack ?
19#define USE_GAUSPIV
20#define USE_LAPACK
21
22// --- Choix de ce que doit faire Lapack
23#define ALSO_LAPACK_INV
24#define ALSO_LAPACK_INV_SYM
25#define ALSO_LAPACK_INV_LSS
26#define ALSO_LAPACK_EV
27#define ALSO_LAPACK_EV_SYM
28
29//////////////////////////////////////////////////
30//////////////////////////////////////////////////
31// NE RIEN CHANGER CI-APRES //
32//////////////////////////////////////////////////
33//////////////////////////////////////////////////
34
35//////////////////////////////////////////////////
[934]36#include "machdefs.h"
[2322]37#include <iostream>
[934]38#include <stdlib.h>
39#include <stdio.h>
40#include <string.h>
41#include <math.h>
42#include <unistd.h>
[1008]43#include "timing.h"
[934]44#include "ntoolsinit.h"
45#include "pexceptions.h"
46#include "array.h"
47#include "srandgen.h"
[2555]48#if defined(USE_LAPACK)
[995]49#include "intflapack.h"
50#endif
[934]51
[2555]52//////////////////////////////////////////////////
[934]53#if defined(COMPLEX)
[975]54 #if defined(PRECIS32)
55 #define TYPE complex<r_4>
[2555]56 #define TYPER r_4
[975]57 #else
58 #define TYPE complex<r_8>
[2555]59 #define TYPER r_8
[975]60 #endif
[2555]61 #define REAL_PART(_x_) (TYPE((_x_).real(),0.))
62 #define CONJ_VAL(_x_) (TYPE((_x_).real(),-(_x_).imag()))
63 #define ABS_VAL(_x_) sqrt((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
[934]64#else
[975]65 #if defined(PRECIS32)
66 #define TYPE r_4
[2555]67 #define TYPER r_4
[975]68 #else
69 #define TYPE r_8
[2555]70 #define TYPER r_8
[975]71 #endif
[2555]72 #define REAL_PART(_x_) (_x_)
73 #define CONJ_VAL(_x_) (_x_)
74 #define ABS_VAL(_x_) fabs((double)_x_)
[934]75#endif
76
[2555]77//////////////////////////////////////////////////
78void Symetrize(TMatrix< TYPE >& A);
79void Hermitian(TMatrix< TYPE >& A);
80r_8 Check_Mat_Ident(TMatrix< TYPE >& A);
81r_8 Check_Mat_VecCol_0(TMatrix< TYPE >& A);
82void Check_Mat_VecCol_2(TMatrix< complex<TYPER> >& A);
83#if defined(USE_LAPACK)
84/*
85-- Pour faire ce test il faut passer la methode ilaenv_en_C()
86 de LapackServer en methode "public" (dans intflapack.h)
87 et recompiler la librairie externe sophya
88*/
89// void TestIlaEnv(int_4 n);
90#endif
91
92
93//////////////////////////////////////////////////
[934]94int main(int narg,char *arg[])
95{
96//--------------------------------------------------------
[2555]97//-- Initialisation
98//--------------------------------------------------------
[934]99// number of lines/columns
[975]100uint_4 N = 5;
[934]101// scale of the value (if =1 values between -1 and 1)
102r_8 scale = 1.;
103// number of values change by +/- vbig
[975]104uint_4 nbig = N;
[2555]105r_8 vbig = 100.;
[975]106// Nombre de lignes de matrice a imprimer
[2555]107uint_4 nprline = N;
[975]108// Initialisation du pauvre de l'aleatoire
[1008]109uint_4 nalea = 0;
110// data scaling for gauss pivoting and determinant
111int tscal = 1;
112bool detok=false;
[2555]113// Please symetrize the input matrice
114bool symetok=false;
115// Please symetrize the input matrice
116bool gaussok=false;
117
[934]118//--------------------------------------------------------
119//-- Decodage arguments
[2555]120//--------------------------------------------------------
[934]121char c;
[2555]122while((c = getopt(narg,arg,"Sdgn:s:b:l:a:t:h")) != -1) {
[934]123 switch (c) {
[2555]124 case 'S' :
125 symetok = true;
[934]126 break;
[2555]127 case 'd' :
128 detok = true;
129 break;
130 case 'g' :
131 gaussok = true;
132 break;
133 case 'n' :
134 sscanf(optarg,"%d",&N);
135 break;
136 case 's' :
137 sscanf(optarg,"%lf",&scale);
138 break;
139 case 'b' :
140 sscanf(optarg,"%d,%lf",&nbig,&vbig);
141 break;
[934]142 case 'l' :
143 sscanf(optarg,"%d",&nprline);
144 break;
[975]145 case 'a' :
146 sscanf(optarg,"%d",&nalea);
147 break;
[2555]148 case 't' :
[1008]149 sscanf(optarg,"%d",&tscal);
150 break;
[934]151 case 'h' :
[2555]152 cout<<"tsttminv [-h] [-n N] [-S] [-s scale] [-b nbig,vbig] [-g]"<<endl
153 <<" [-l nprline] [-a nalea] [-t tscal] [-d]"<<endl
154 <<"-- matrix A(N,N) filled with {[-1,1] +/- vbig(nbig time)}*scale --"<<endl
155 <<"-g : instead of flat [-1,1] use normal gaussian distribution for A(i,j)"<<endl
156 <<"-S : symetrize the input matrix"<<endl
157 <<"-l : print nprline of input and test matrices"<<endl
158 <<"-a : for random (pseudo) changing"<<endl
159 <<"-- Only GausPiv --"<<endl
160 <<"-t 0/1/2 : data scaling 0=no, 1=global (def), 2=row-by-row"<<endl
[2557]161 <<"-d : also compute determinant (ne pas utiliser si N est grand)"<<endl;
[2555]162 return(-1);
[934]163 }
164}
165if(N<=1) N = 1;
[975]166cout<<"Taille matrice NxN, N = "<<N<<endl;
[2555]167if(gaussok) cout<<"Elements gaussian normal * scale = "<<scale<<endl;
168 else cout<<"Elements entre +/- 1 * scale = "<<scale<<endl;
[975]169cout<<"Nombre de valeurs hors standard nbig = "<<nbig<<endl;
[2555]170cout<<"Valeurs hors standard (+/- vbig = "<<vbig<<" ) * scale = "<<vbig*scale<<endl;
[934]171cout<<"Nombre de lignes de matrice a imprimer "<<nprline<<endl;
[975]172cout<<"Initialisation de l aleatoire par "<<nalea<<" tirages"<<endl;
[2555]173cout<<"Data scaling "<<tscal<<" determinant="<<detok<<endl;
174if(symetok) cout<<"Input matrix has been symetrized "<<symetok<<endl;
[934]175cout<<endl;
176
[2555]177//--------------------------------------------------------
178// TestIlaEnv(N); return -41;
179//--------------------------------------------------------
180
181//--------------------------------------------------------
[934]182//-- Initialization arrays
[2555]183//--------------------------------------------------------
[934]184SophyaInit();
[2555]185InitTim();
186#if defined(USE_LAPACK)
[995]187 BaseArray::SetDefaultMemoryMapping(BaseArray::FortranMemoryMapping);
188#endif
[2555]189if(nalea>0) for(int i=0;i<nalea;i++) drand01();
190BaseArray::SetMaxPrint(nprline*N,0);
[934]191
[2555]192//--------------------------------------------------------
193//-- Definition global arrays
194//--------------------------------------------------------
195TMatrix< TYPE > Ainput(N,N); Ainput = (TYPE) 0;
196TMatrix< TYPE > A(N,N); A = (TYPE) 0;
197Ainput.Show();
[934]198
[2555]199//--------------------------------------------------------
200//-- Fill matrices with flat random
201//--------------------------------------------------------
202if(gaussok) Ainput = RandomSequence(RandomSequence::Gaussian,0.,1.);
203 else Ainput = RandomSequence(RandomSequence::Flat,0.,1.);
204#if defined(COMPLEX)
205if(gaussok) A = RandomSequence(RandomSequence::Gaussian,0.,1.);
206 else A = RandomSequence(RandomSequence::Flat,0.,1.);
207Ainput += TYPE(0.,1.)*A;
208#endif
[1008]209
[2555]210//--------------------------------------------------------
211//-- Fill matrices with big values
212//--------------------------------------------------------
213if(nbig>0) {
214#if defined(COMPLEX)
215 nbig = (nbig+1)/2;
216#endif
217 TMatrix< uint_2 > Vind(N,N); Vind = 0;
218 // for real part
219 uint_4 nbr=0;
220 for(int k=0;k<nbig;k++) {
221 int i = (int) (drand01()*N); int j = (int) (drand01()*N);
222 double s=(drand01()>0.5)?1.:-1.;
223 if(Vind(i,j)==0) {Ainput(i,j) += (TYPER) s*vbig; Vind(i,j)+=1; nbr++;}
224 }
225 cout<<"Nombre de valeurs BIG reelles = "<<nbr<<endl;
226#if defined(COMPLEX)
227 // for imaginary part
228 uint_4 nbi=0;
229 for(int k=0;k<nbig;k++) {
230 int i = (int) (drand01()*N); int j = (int) (drand01()*N);
231 double s=(drand01()>0.5)?1.:-1.;
232 if(Vind(i,j)<=1) {Ainput(i,j) += TYPE(0.,(TYPER)s*vbig); Vind(i,j)+=2; nbi++;}
233 }
234 cout<<"Nombre de valeurs BIG imaginaires = "<<nbi<<endl;
235 cout<<"Nombre de valeurs BIG = "<<nbr+nbi<<endl;
236#endif
237}
[934]238
[2555]239//--------------------------------------------------------
240//-- Scale matrix for machine precision tests
241//--------------------------------------------------------
242Ainput *= (TYPE) scale;
243
244//--------------------------------------------------------
245//-- Create symetric matrix for all A if requested
246//--------------------------------------------------------
247if(symetok) Symetrize(Ainput);
248
249//--------------------------------------------------------
250//-- Print matrice Ainput
251//--------------------------------------------------------
252cout<<"------------ TMatrix Ainput :"<<endl;
253if(nprline>0) {cout<<Ainput; cout<<endl;}
254PrtTim("--- End of Matrix filling ---");
255
256
[2557]257//////////////////////////////////////////////
258///////// Test Inversion avec Lapack /////////
259//////////////////////////////////////////////
260#if defined(USE_LAPACK) && defined(ALSO_LAPACK_INV)
[2555]261{
262cout<<"\n=========================================="<<endl;
263cout<<"------------ Inversion LAPACK"<<endl;
264A = Ainput;
265//-- Inversion
266TMatrix< TYPE > InvA(N,N); InvA = IdentityMatrix(1.,N);
267int_4 info = LapackLinSolve(A,InvA);
268cout<<"info="<<info<<endl;
269PrtTim("--- End of LapackLinSolve Inversion ---");
270//-- AiA = A * InvA
271cout<<"Compute AiA = A * InvA"<<endl;
272TMatrix< TYPE > AiA(N,N); AiA = Ainput * InvA;
273cout<<"------------ TMatrix AiA = A * InvA:"<<endl;
274if(nprline>0) {cout<<AiA; cout<<endl;}
275//-- Check
276Check_Mat_Ident(AiA);
277PrtTim("--- End of LapackLinSolve Test ---");
[975]278}
[2555]279#endif
280
281
[2557]282//////////////////////////////////////////////////
283///////// Test Inversion avec Lapack sym /////////
284//////////////////////////////////////////////////
285#if defined(USE_LAPACK) && defined(ALSO_LAPACK_INV_SYM)
[2555]286{
287cout<<"\n=========================================="<<endl;
288cout<<"------------ Inversion LAPACK sym"<<endl;
289TMatrix< TYPE > Asym(N,N); Asym=Ainput; Symetrize(Asym); A=Asym;
290//-- Inversion
291TMatrix< TYPE > InvA(N,N); InvA = IdentityMatrix(1.,N);
292int_4 info = LapackLinSolveSym(A,InvA);
293cout<<"info="<<info<<endl;
294PrtTim("--- End of LapackLinSolveSym Inversion ---");
295//-- AiA = A * InvA
296cout<<"Compute AiA = A * InvA"<<endl;
297TMatrix< TYPE > AiA(N,N); AiA = Asym * InvA;
298cout<<"------------ TMatrix AiA = A * InvA:"<<endl;
299if(nprline>0) {cout<<AiA; cout<<endl;}
300//-- Check
301Check_Mat_Ident(AiA);
302PrtTim("--- End of LapackLinSolveSym Test ---");
[975]303}
[934]304#endif
305
306
[2555]307////////////////////////////////////////////////
308///////// Test avec Lapack LeastSquare /////////
309////////////////////////////////////////////////
[2557]310#if defined(USE_LAPACK) && defined(ALSO_LAPACK_INV_LSS)
[2555]311{
312cout<<"\n=========================================="<<endl;
313cout<<"------------ Inversion LAPACK LeastSquare"<<endl;
314A = Ainput;
[934]315//-- Inversion
[2555]316TMatrix< TYPE > InvA(N,N); InvA = IdentityMatrix(1.,N);
317int_4 info = LapackLeastSquareSolve(A,InvA);
318cout<<"info="<<info<<endl;
319PrtTim("--- End of LapackLeastSquareSolve Inversion ---");
320//-- AiA = A * InvA
321cout<<"Compute AiA = A * InvA"<<endl;
322TMatrix< TYPE > AiA(N,N); AiA = Ainput * InvA;
323cout<<"------------ TMatrix AiA = A * InvA:"<<endl;
324if(nprline>0) {cout<<AiA; cout<<endl;}
325//-- Check
326Check_Mat_Ident(AiA);
327PrtTim("--- End of LapackLeastSquareSolve Test ---");
328}
329#endif
330
331
332////////////////////////////////////////////////
333///////// Test avec Lapack pour EV /////////
334////////////////////////////////////////////////
[2557]335#if defined(USE_LAPACK) && defined(ALSO_LAPACK_EV)
[2555]336{
337cout<<"\n=========================================="<<endl;
338cout<<"------------ Eigen decompositon LapackEigen "<<endl;
339A=Ainput;
340TMatrix< TYPE > Evec(N,N); Evec = (TYPE) 0;
341TVector< complex<r_8> > Eval(N); Eval = complex<r_8>(0,0);
342//-- Decompositon
343int_4 info = LapackEigen(A,Eval,Evec,true);
344cout<<"info="<<info<<endl;
345PrtTim("--- End of LapackEigen decompositon ---");
346if(nprline>0) {cout<<Eval; cout<<endl; cout<<Evec; cout<<endl;}
347//-- Find the complex conjugate pairs
[2557]348#if !defined(COMPLEX)
[2555]349TVector< uint_2 > Evalconj(N); Evalconj = 0;
350int_4 nconj=0;
351for(int i=0;i<N-1;i++) {
352 if(Evalconj(i)!=0) continue; // deja traite
353 if(Eval(i).imag()==0.) continue; // real eigenvalue
354 if(fabs(Eval(i).imag()+Eval(i+1).imag())>1e-150) continue; // les 2 consecutives ne sont pas conjuguees
355 if(Eval(i).imag()<0.) continue; // first conjugate have positive imaginary part
356 if(Eval(i+1).imag()>0.) continue;
357 Evalconj(i) = 1; Evalconj(i+1) = 2; nconj++;
358}
359//cout<<Evalconj<<endl;
360cout<<"Number of conjugate eigen values: "<<nconj<<" *2 = "<<2*nconj<<" / "<<N<<endl;
361#endif
362//-- Azmlz = A*z(l) - l*z(l)
363cout<<"Compute Azmlz(l) = A*z(l) - l*z(l)"<<endl;
364TMatrix< complex<TYPER> > Azmlz(N,N); Azmlz = (complex<TYPER>) 0;
365for(int l=0;l<N;l++) { // eigen value
366 complex<TYPER> Eval_l = complex<TYPER>(Eval(l).real(),Eval(l).imag());
367 for(int i=0;i<N;i++) {
368 complex<TYPER> Evec_il;
[2557]369 #if defined(COMPLEX)
[2555]370 Evec_il = Evec(i,l);
371 #else
372 Evec_il = complex<TYPER>(Evec(i,l),0.);
373 if(Evalconj(l)==1) Evec_il = complex<TYPER>(Evec(i,l),Evec(i,l+1));
374 else if(Evalconj(l)==2) Evec_il = complex<TYPER>(Evec(i,l-1),-Evec(i,l));
375 #endif
376 for(int j=0;j<N;j++) {
377 complex<TYPER> Evec_jl;
[2557]378 #if defined(COMPLEX)
[2555]379 Evec_jl = Evec(j,l);
380 #else
381 Evec_jl = complex<TYPER>(Evec(j,l),0.);
382 if(Evalconj(l)==1) Evec_jl = complex<TYPER>(Evec(j,l),Evec(j,l+1));
383 else if(Evalconj(l)==2) Evec_jl = complex<TYPER>(Evec(j,l-1),-Evec(j,l));
384 #endif
385 Azmlz(i,l) += Ainput(i,j) * Evec_jl;
386 }
387 Azmlz(i,l) -= Eval_l*Evec_il;
388 }
389}
390if(nprline>0) {cout<<Azmlz; cout<<endl;}
391//-- Check
392Check_Mat_VecCol_2(Azmlz);
393PrtTim("--- End of LapackEigen Test ---");
394}
395#endif
396
397
398////////////////////////////////////////////////
399///////// Test avec Lapack sym pour EV /////////
400////////////////////////////////////////////////
[2557]401#if defined(USE_LAPACK) && defined(ALSO_LAPACK_EV_SYM)
[2555]402{
403cout<<"\n=========================================="<<endl;
404cout<<"------------ Eigen decompositon LapackEigenSym "<<endl;
405TMatrix< TYPE > Aher(N,N); Aher=Ainput; Hermitian(Aher); A=Aher;
406TVector<r_8> Eval;
407//-- Decompositon
408int_4 info = LapackEigenSym(A,Eval,true);
409cout<<"info="<<info<<endl;
410PrtTim("--- End of LapackEigenSym decompositon ---");
411if(nprline>0) {cout<<Eval; cout<<endl; cout<<A; cout<<endl;}
412//-- Azmlz = A*z(l) - l*z(l)
413// le vecteur propre z pour la l-ieme valeur propre est dans A(.,l):
414// z_i = A(i,l) ou "l" est la l-ieme valeur propre
415cout<<"Compute Azmlz(l) = A*z(l) - l*z(l)"<<endl;
416TMatrix< TYPE > Azmlz(N,N); Azmlz = (TYPE) 0;
417for(int l=0;l<N;l++) // eigen value
418 for(int i=0;i<N;i++)
419 {for(int j=0;j<N;j++) Azmlz(i,l)+=Aher(i,j)*A(j,l); Azmlz(i,l)-=(TYPER)Eval(l)*A(i,l);}
420if(nprline>0) {cout<<Azmlz; cout<<endl;}
421//-- Check
422Check_Mat_VecCol_0(Azmlz);
423PrtTim("--- End of LapackEigenSym Test ---");
424}
425#endif
426
427
[2557]428///////////////////////////////////////////////
429///////// Test Inversion avec GausPiv /////////
430///////////////////////////////////////////////
431#if defined(USE_GAUSPIV)
[2555]432{
433cout<<"\n==========================================\n"
434 <<"------------ Inversion GausPiv"<<endl;
435SimpleMatrixOperation< TYPE >::SetGausPivScal(tscal);
436A = Ainput;
437//-- Inversion
438TMatrix< TYPE > InvA(N,N); InvA = IdentityMatrix(1.,N);
[1008]439TYPE det = GausPiv(A,InvA,detok);
[2555]440PrtTim("--- End of GausPiv Inversion ---");
[1008]441cout<<"Det = "<<det<<endl;
[934]442cout<<"------------ TMatrix InvA = A^(-1):"<<endl;
443//-- AiA = A * InvA
[2555]444cout<<"Compute AiA = A * InvA"<<endl;
445TMatrix< TYPE > AiA(N,N); AiA = Ainput * InvA;
[934]446cout<<"------------ TMatrix AiA = A * InvA:"<<endl;
447if(nprline>0) {cout<<AiA; cout<<endl;}
448//-- Check
[2555]449Check_Mat_Ident(AiA);
450PrtTim("--- End of GausPiv Test ---");
[934]451}
[2555]452#endif
453
454
455PrtTim("--- End of Job ---");
456exit(0);
[934]457}
458
[995]459
[2555]460
461////////////////////////////////////////////////////////////
462////-------------------------------------------------------
463void Symetrize(TMatrix< TYPE >& A)
464// Symetrize A
465{
466 int_4 N = A.NRows();
467 for(int i=0;i<N-1;i++) for(int j=i+1;j<N;j++) A(j,i) = A(i,j);
[995]468}
[2555]469
470////-------------------------------------------------------
471void Hermitian(TMatrix< TYPE >& A)
472// Put A hermitian
473{
474 int_4 N = A.NRows();
475 for(int i=0;i<N-1;i++) for(int j=i+1;j<N;j++) A(j,i) = CONJ_VAL(A(i,j));
476 for(int i=0;i<N;i++) A(i,i) = REAL_PART(A(i,i));
[995]477}
478
[2555]479////-------------------------------------------------------
480r_8 Check_Mat_Ident(TMatrix< TYPE >& A)
481// Compute the biggest difference element by element of A / Identity
482{
483 int_4 N = A.NRows();
484 r_8 vmaxd=-1.;
485 for(int i=0;i<N;i++)
486 if( ABS_VAL((TYPER)1.-A(i,i)) > vmaxd ) vmaxd = ABS_VAL((TYPER)1.-A(i,i));
487 cout<<"Ecart maximum par rapport a 1 sur diagonale = "<<vmaxd<<endl;
488 r_8 vmaxh = -1.;
489 for(int i=0;i<N;i++) for(int j=0;j<N;j++) {
490 if(i==j) continue;
491 if( ABS_VAL(A(i,j)) > vmaxh ) vmaxh = ABS_VAL(A(i,j));
492 }
493 cout<<"Ecart maximum par rapport a 0 hors diagonale = "<<vmaxh<<endl;
494 return (vmaxd>vmaxh)? vmaxd: vmaxh;
495}
[995]496
[2555]497////-------------------------------------------------------
498r_8 Check_Mat_VecCol_0(TMatrix< TYPE >& A)
499// Return the biggest norm of the N vectors column of matrix
500{
501 int_4 N = A.NRows();
502 r_8 vmax=-1.;
503 for(int l=0;l<N;l++) {
504 r_8 absv = 0.;
505 for(int i=0;i<N;i++) absv += ABS_VAL(A(i,l)) * ABS_VAL(A(i,l));
506 if( absv > vmax ) vmax = absv;
507 }
508 vmax = sqrt(vmax);
509 cout<<"Longueur max de ||A*z-l*z|| pour tous l = "<<vmax<<endl;
510 return vmax;
[934]511}
[2555]512
513////-------------------------------------------------------
514void Check_Mat_VecCol_2(TMatrix< complex<TYPER> >& A)
515// Return the biggest norm of :
516// - the real part of the N vectors column of matrix
517// - the imaginary part of the N vectors column of matrix
518// - the module of the N vectors column of matrix
519{
520 int_4 N = A.NRows();
521 r_8 vmaxr=-1., vmaxi=-1., vmaxn=-1.;
522 for(int l=0;l<N;l++) {
523 double absvr = 0., absvi = 0., absvn = 0.;
524 for(int i=0;i<N;i++) {
525 absvr += A(i,l).real()*A(i,l).real();
526 absvi += A(i,l).imag()*A(i,l).imag();
527 absvn += A(i,l).real()*A(i,l).real() + A(i,l).imag()*A(i,l).imag();
528 }
529 if( absvr > vmaxr ) vmaxr = absvr;
530 if( absvi > vmaxi ) vmaxi = absvi;
531 if( absvn > vmaxn ) vmaxn = absvn;
532 }
533 vmaxr=sqrt(vmaxr); vmaxi=sqrt(vmaxi); vmaxn=sqrt(vmaxn);
534 cout<<"Longueur max de ||A*z-l*z|| pour tous l, reel = "<<vmaxr
535 <<", imag = "<<vmaxi<<", module = "<<vmaxn<<endl;
536}
537
538
539/*
540void TestIlaEnv(int_4 n)
541{
542LapackServer<TYPE> lps;
543cout<<"TestIlaEnv n="<<n<<endl;
544cout<<lps.ilaenv_en_C(1,"SSYTRF","U",n,-1,-1,-1)<<endl;
545cout<<lps.ilaenv_en_C(1,"DSYTRF","U",n,-1,-1,-1)<<endl;
546cout<<lps.ilaenv_en_C(1,"CSYTRF","U",n,-1,-1,-1)<<endl;
547cout<<lps.ilaenv_en_C(1,"ZSYTRF","U",n,-1,-1,-1)<<endl;
548cout<<lps.ilaenv_en_C(1,"SSYTRF","L",n,-1,-1,-1)<<endl;
549cout<<lps.ilaenv_en_C(1,"DSYTRF","L",n,-1,-1,-1)<<endl;
550cout<<lps.ilaenv_en_C(1,"CSYTRF","L",n,-1,-1,-1)<<endl;
551cout<<lps.ilaenv_en_C(1,"ZSYTRF","L",n,-1,-1,-1)<<endl;
552}
553*/
Note: See TracBrowser for help on using the repository browser.