-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmain.c
930 lines (771 loc) · 38.1 KB
/
main.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
/**
*
* Copyright (c) 2017, King Abdullah University of Science and Technology
* All rights reserved.
*
**/
/**
*
* @file testing_pdgeqdwh.c
*
* QDWH is a high performance software framework for computing
* the polar decomposition on distributed-memory manycore systems provided by KAUST
*
* @version 2.0.0
* @author Dalal Sukkari
* @author Hatem Ltaief
* @date 2017-11-13
*
**/
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <assert.h>
#include <math.h>
#include <getopt.h>
#include <unistd.h>
#include <sys/time.h>
#include <sys/resource.h>
#include <mpi.h>
#include <mkl_lapack.h>
#include <mkl_lapacke.h>
#include "myscalapack.h"
#include "include.h"
#include "flops.h"
/* Default values of parameters */
int nprow = 1;
int npcol = 1;
int lvec = 1;
int rvec = 1;
int m = 5120;
int n = 5120;
int nb = 128;
int mode = 0;
double cond = 9.0072e+15;
int optcond = 0;
int start = 5120;
int stop = 5120;
int step = 1;
int niter = 1;
int check = 0;
int verbose = 0;
void print_usage(void)
{
fprintf(stderr,
"======= QDWH testing using ScaLAPACK\n"
" -p --nprow : Number of MPI process rows\n"
" -q --npcol : Number of MPI process cols\n"
" -jl --lvec : Compute left singular vectors\n"
" -jr --rvec : Compute right singular vectors\n"
" -m --M : First dimension of the matrix\n"
" -n --N : Second dimension of the matrix\n"
" -b --nb : Block size\n"
" -d --mode : [1:6] Mode from pdlatms used to generate the matrix\n"
" -k --cond : Condition number used to generate the matrix\n"
" -o --optcond : Estimate Condition number using QR\n"
" -i --niter : Number of iterations\n"
" -r --n_range : Range for matrix sizes Start:Stop:Step\n"
" -c --check : Check the solution\n"
" -v --verbose : Verbose\n"
" -h --help : Print this help\n" );
}
#define GETOPT_STRING "p:q:x:y:m:n:b:d:i:o:r:Q,S:s:w:e:c:f:t:v:h"
static struct option long_options[] =
{
/* PaRSEC specific options */
{"nprow", required_argument, 0, 'p'},
{"npcol", required_argument, 0, 'q'},
{"jl", no_argument, 0, 'x'},
{"lvec", no_argument, 0, 'x'},
{"jr", no_argument, 0, 'y'},
{"rvec", no_argument, 0, 'y'},
{"M", required_argument, 0, 'm'},
{"m", required_argument, 0, 'm'},
{"N", required_argument, 0, 'n'},
{"n", required_argument, 0, 'n'},
{"nb", required_argument, 0, 'b'},
{"b", required_argument, 0, 'b'},
{"mode", required_argument, 0, 'd'},
{"d", required_argument, 0, 'd'},
{"cond", required_argument, 0, 'k'},
{"k", required_argument, 0, 'k'},
{"optcond", required_argument, 0, 'o'},
{"o", required_argument, 0, 'o'},
{"i", required_argument, 0, 'i'},
{"niter", required_argument, 0, 'i'},
{"r", required_argument, 0, 'r'},
{"n_range", required_argument, 0, 'r'},
{"check", no_argument, 0, 'c'},
{"verbose", no_argument, 0, 'v'},
{"help", no_argument, 0, 'h'},
{"h", no_argument, 0, 'h'},
{0, 0, 0, 0}
};
static void parse_arguments(int argc, char** argv)
{
int opt = 0;
int c;
int myrank_mpi;
MPI_Comm_rank(MPI_COMM_WORLD, &myrank_mpi);
do {
#if defined(HAVE_GETOPT_LONG)
c = getopt_long_only(argc, argv, "",
long_options, &opt);
#else
c = getopt(argc, argv, GETOPT_STRING);
(void) opt;
#endif /* defined(HAVE_GETOPT_LONG) */
switch(c) {
case 'p': nprow = atoi(optarg); break;
case 'q': npcol = atoi(optarg); break;
case 'm': n = atoi(optarg); start = m; stop = m; step = 1; break;
case 'n': n = atoi(optarg); //start = n; stop = n; step = 1; break;
case 'b': nb = atoi(optarg); break;
case 'd': mode = atoi(optarg); break;
case 'k': cond = atof(optarg); break;
case 'o': optcond = atof(optarg); break;
case 'i': niter = atoi(optarg); break;
case 'r': get_range( optarg, &start, &stop, &step ); break;
case 'c': check = 1; break;
case 'v': verbose = 1; break;
case 'h':
if (myrank_mpi == 0) print_usage(); MPI_Finalize(); exit(0);
break;
default:
break;
}
} while(-1 != c);
}
int main(int argc, char **argv) {
int myrank_mpi, nprocs_mpi;
int ictxt, myrow, mycol;
int mloc, nloc, mlocW;
int mloc_min_mn, nloc_min_mn;
int mloc_max_mn, nloc_max_mn;
int mloc_n, nloc_n;
int mloc_pinv, nloc_pinv;
int mpi_comm_rows, mpi_comm_cols;
int i, j, iter, size, info_facto_mr, info_facto_dc, info_facto_qw, info_facto_el, info_facto_sl, info, iseed;
int my_info_facto;
int i0 = 0, i1 = 1;
int lwork, liwork, ldw, lwork_pdgesvd;
/* Quick return if possible */
if ( m == 0 || n == 0 ) {
return 0;
}
/* Check the inputs */
if ( m < 0 || n < 0 ){
printf("\n matrix with negative dimension \n");
return 1;
}
/* Allocation for the input/output matrices */
//int descA[9], descH[9];
//double *A=NULL, *H=NULL;
/* Allocation to check the results */
//int descAcpy[9], descC[9];
//double *Acpy=NULL, *C=NULL;
/* Allocation for pdlatsm */
//double *Wloc1=NULL, *Wloc2=NULL, *D=NULL;
double eps = LAPACKE_dlamch_work('e');
int iprepad, ipostpad, sizemqrleft, sizemqrright, sizeqrf, sizeqtq,
sizechk, sizesyevx, isizesyevx,
sizesubtst, isizesubtst, sizetst,
isizetst;
double *D, *Wloc1;
/**/
double flops, GFLOPS;
double orth_Uqw, berr_UHqw;
double frobA;
double alpha, beta;
char *jobu, *jobvt;
/**/
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Program starts... \n");
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &myrank_mpi);
MPI_Comm_size(MPI_COMM_WORLD, &nprocs_mpi);
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "MPI Init done\n");
parse_arguments(argc, argv);
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Checking arguments done\n");
Cblacs_get( -1, 0, &ictxt );
Cblacs_gridinit( &ictxt, "R", nprow, npcol );
Cblacs_gridinfo( ictxt, &nprow, &npcol, &myrow, &mycol );
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "BLACS Init done\n");
if ( myrank_mpi == 0 ) {
fprintf(stderr, "# \n");
fprintf(stderr, "# NPROCS %d P %d Q %d\n", nprocs_mpi, nprow, npcol);
fprintf(stderr, "# niter %d\n", niter);
fprintf(stderr, "# n_range %d:%d:%d mode: %d cond: %2.4e \n", start, stop, step, mode, cond);
fprintf(stderr, "# \n");
}
//int i, j;
//int info;
double my_elapsed_qwsvd, elapsed_qwsvd, sumtime_qwsvd, min_time_qwsvd, max_time_qwsvd;
double my_elapsed_pdgesvd, elapsed_pdgesvd, sumtime_pdgesvd, min_time_pdgesvd, max_time_pdgesvd, flops_pdgesvd;
long int min_mn, max_mn, n2;
//int n232 = r32up(n2);
//int lddb = n232;
//int ldd_min_mn = r32up(min_mn);
//int ldd_max_mn = r32up(max_mn);
//int ldda = ldd_max_mn;
//int m32 = r32up(m);
//int n32 = r32up(n);
//int lddu = m32;
//int lddvt = n32;
/*
* Allocatation for PDGESVD
*/
double *Work_pdgesvd;
/*
* Allocatation for QDWHPartial
*/
//float *d_A, *d_B, *d_VT, *d_U, *d_pinv, *S;
int descA[9], descB[9], descVT[9], descU[9], descpinv[9], descS[9];
double *A=NULL, *B=NULL, *VT=NULL, *U=NULL, *pinv=NULL, *S=NULL;
/*
* Allocatation to check the results
*/
int descAcpy[9];
double *Acpy=NULL, *W=NULL;
int sizeS = 0, k = 0, ktmp, sizeStmp;
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Range loop starts\n");
// Begin loop over range
for (size = start; size <= stop; size += step) {
while ( (int)((double)size / (double)nb) < ( max(nprow , npcol) )){
if ( myrank_mpi == 0 ) fprintf(stderr, " Matrix size is small to be facrorized using this number of processors \n");
size += step;
}
m = size;
n = size;
ldw = 2*min(m,n);
min_mn = min(m,n);
max_mn = max(m,n);
n2 = 2*min_mn;
my_elapsed_qwsvd = 0.0, elapsed_qwsvd = 0.0, sumtime_qwsvd = 0.0, min_time_qwsvd = 0.0, max_time_qwsvd = 0.0;
my_elapsed_pdgesvd = 0.0, elapsed_pdgesvd = 0.0, sumtime_pdgesvd = 0.0, min_time_pdgesvd = 0.0, max_time_pdgesvd = 0.0;
mloc = numroc_( &m, &nb, &myrow, &i0, &nprow );
nloc = numroc_( &n, &nb, &mycol, &i0, &npcol );
mloc_min_mn = numroc_( &min_mn, &nb, &myrow, &i0, &nprow );
nloc_min_mn = numroc_( &min_mn, &nb, &mycol, &i0, &npcol );
mloc_max_mn = numroc_( &max_mn, &nb, &myrow, &i0, &nprow );
nloc_max_mn = numroc_( &max_mn, &nb, &mycol, &i0, &npcol );
mloc_n = numroc_( &n, &nb, &myrow, &i0, &nprow );
nloc_n = numroc_( &n, &nb, &mycol, &i0, &npcol );
mloc_pinv = numroc_( &n, &nb, &myrow, &i0, &nprow );
nloc_pinv = numroc_( &m, &nb, &mycol, &i0, &npcol );
mlocW = numroc_( &ldw, &nb, &myrow, &i0, &nprow );
//printf("\n mloc %d nloc %d mloc_min_mn %d nloc_min_mn %d mloc_max_mn %d nloc_max_mn %d mlocW %d mloc_pinv %d nloc_pinv %d \n", mloc, nloc, mloc_min_mn, nloc_min_mn, mloc_max_mn, nloc_max_mn, mlocW, mloc_pinv, nloc_pinv);
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Desc Init starts %d\n", mloc);
/*
* Allocatation for QDWHPartial
*/
sizeS = 0, k = 0;
descinit_( descA, &m, &n, &nb, &nb, &i0, &i0, &ictxt, &mloc, &info );
descinit_( descB, &ldw, &min_mn, &nb, &nb, &i0, &i0, &ictxt, &mlocW, &info );
descinit_( descU, &m, &n, &nb, &nb, &i0, &i0, &ictxt, &mloc, &info );
//descinit_( descU, &m, &min_mn, &nb, &nb, &i0, &i0, &ictxt, &mloc, &info );
descinit_( descVT, &n, &min_mn, &nb, &nb, &i0, &i0, &ictxt, &mloc_n, &info );
descinit_( descpinv, &n, &m, &nb, &nb, &i0, &i0, &ictxt, &mloc_pinv, &info );
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Desc Init ends %d\n", mloc);
A = (double *)malloc(mloc*nloc*sizeof(double));
B = (double *)malloc(mlocW*nloc_min_mn*sizeof(double));
//U = (double *)malloc(mloc*nloc_min_mn*sizeof(double)) ;
U = (double *)malloc(mloc*nloc*sizeof(double));
VT = (double *)malloc(mloc_n*nloc_min_mn*sizeof(double));
pinv = (double *)malloc(mloc_pinv*nloc_pinv*sizeof(double));
S = (double *)malloc(min_mn*sizeof(double));
/*
* Allocatation to check the results
float *A, *Acpy, *test_pinv1, *test_pinv2, *pinv;
magma_smalloc_pinned( &A, m*n );
magma_smalloc_pinned( &Acpy, m*n );
magma_smalloc_pinned( &test_pinv1, m*m );
magma_smalloc_pinned( &test_pinv2, m*n );
magma_smalloc_pinned( &pinv, n*m );
lapackf77_slaset( "A", &m, &n, &alpha, &alpha, A, &m );
lapackf77_slaset( "A", &m, &n, &alpha, &alpha, Acpy, &m );
*/
Acpy = (double *)malloc(mloc*nloc*sizeof(double)) ;
descinit_( descAcpy, &m, &n, &nb, &nb, &i0, &i0, &ictxt, &mloc, &info );
/*
* Allocatation for DLATMS
*/
D = (double *)malloc(min_mn*sizeof(double)) ;
W = (double *)malloc(min_mn*sizeof(double)) ;
char *dist = "N"; /* NORMAL( 0, 1 ) ( 'N' for normal ) */
int iseed[4] = {1, 0, 0, 1};
char *sym = "P"; /* The generated matrix is symmetric, with
eigenvalues (= singular values) specified by D, COND,
MODE, and DMAX; they will not be negative.
"N" not supported. */
//int mode = 4; /* sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND) */
//double cond = 1.0/eps;
double dmax = 1.0;
int kl = n;
int ku = n;
char *pack = "N"; /* no packing */
int order = n;
int info;
pdlasizesep_( descA,
&iprepad, &ipostpad, &sizemqrleft, &sizemqrright, &sizeqrf,
&lwork,
&sizeqtq, &sizechk, &sizesyevx, &isizesyevx, &sizesubtst,
&isizesubtst, &sizetst, &isizetst );
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Setting lwork done\n");
Wloc1 = (double *)calloc(max(m,lwork),sizeof(double)) ;
/*
* For matrix similar to matlab
* Generate random synthetic matrices using SLATMS
* Condition number = sigma_max / Sigma_min (if the mode = 0)
*/
mode = 0;
//double cond = 1.0/eps*1.e0;
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Setting singular values \n");
if ( myrank_mpi == 0 ) printf("\n min_mn %d \n", min_mn);
for( i = 0; i < min_mn; i++ ){
W[i] = pow(0.5, (double)i/(double)min_mn*100.); // geometrically increasing
}
double cond = W[0] / W[min_mn-1];
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Setting singular values done\n");
pdlatms_(&m, &n, dist,
iseed, sym, W, &mode, &cond, &dmax,
&m, &n, pack,
A, &i1, &i1, descA, &n,
Wloc1, &lwork, &info);
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "MatGen done\n");
if ( info != 0 ) {
fprintf(stderr, "An error occured during matrix generation: %d\n", info );
return EXIT_FAILURE;
}
/*
* Save copy of the matrix
*/
pdlacpy_( "A", &m, &n,
A, &i1, &i1, descA,
Acpy, &i1, &i1, descAcpy );
{
/*
* Find SVD using SGESVD
double *Wsvd, *Usvd, *Vsvd;
int descUsvd[9], descVsvd[9];
int lworksvd = -1;
Wsvd = (double *)malloc(1*sizeof(double));
pdgesvd_( "N", "N", &m, &n, A, &i1, &i1, descA,
W,
Usvd, &i1, &i1, descUsvd,
Vsvd, &i1, &i1, descVsvd,
Wsvd, &lworksvd, &my_info_facto );
lworksvd = Wsvd[0];
Wsvd = (double *)malloc(lworksvd*sizeof(double));
pdgesvd_( "N", "N", &m, &n, A, &i1, &i1, descA,
W,
Usvd, &i1, &i1, descUsvd,
Vsvd, &i1, &i1, descVsvd,
Wsvd, &lworksvd, &my_info_facto );
free(Wsvd);
*/
}
double normA = pdlange_ ( "f", &m, &n, Acpy, &i1, &i1, descAcpy, Wloc1);
// loop over iteration
for ( iter = 0; iter < niter; iter++ ) {
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "\nScaLAPACK dgesvd starts...\n");
my_elapsed_qwsvd = 0.0, elapsed_qwsvd = 0.0;
my_elapsed_pdgesvd = 0.0, elapsed_pdgesvd = 0.0;
/*
* Save copy of A in Acpy
*/
pdlacpy_( "A", &m, &n,
Acpy, &i1, &i1, descAcpy,
A, &i1, &i1, descA );
/************************************************************
* Call QDWHPartial
************************************************************/
{
double s = 1.e-4; //Threshold for the singular values
double tol = 1.e-2; //Tolerance (governs accuracy)
int it = 0;
int fact = 1;
int psinv = 0;
flops = 0.0;
if ( myrank_mpi == 0 ) printf("\n Start QDWHpartial ... \n");
my_elapsed_qwsvd = 0.0;
my_elapsed_qwsvd =- MPI_Wtime();
QDWHpartial ( m, n, // Size of matrix
fact, psinv,
s, // Threshold
tol, // Tolerance (governs accuracy)
A, i1, i1, descA, // Matrix
S, // Sigular values, size min(m,n)
U, i1, i1, descU, // Left singular vectors, size mx(10%n)
VT, i1, i1, descVT, // Right singular vectors, size nxn, d_VT = d_VT
B, i1, i1, descB,// Needed for the QR fact in QDWH, it is of size NxN, because the matrix will reduced
&sizeS, // Num-computed-SV
&k, // Num-wanted-SV
&it,
&flops);
my_elapsed_qwsvd += MPI_Wtime();
MPI_Allreduce( &my_elapsed_qwsvd, &elapsed_qwsvd, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
sumtime_qwsvd += elapsed_qwsvd;
if ( elapsed_qwsvd >= max_time_qwsvd ) { max_time_qwsvd = elapsed_qwsvd;}
if ( elapsed_qwsvd <= min_time_qwsvd ) { min_time_qwsvd = elapsed_qwsvd;}
if ( myrank_mpi == 0 ) printf("\n Done QDWHpartial ... \n");
ktmp = k;
sizeStmp = sizeS;
k = min(ktmp,sizeStmp);
sizeS = k;
/*
* Check the results
*/
if ( check )
{
if ( myrank_mpi == 0 ) {
fprintf(stderr, "\n QDWHpartial \n");
fprintf(stderr, "/////////////////////////////////////////////////////////////////////////\n");
}
for (i = 0; i < k; i++){
D[i] = W[i] - S[i];
}
alpha = 1.0; beta = -1.0;
int ione = 1;
double acc_sv = LAPACKE_dlange_work( LAPACK_COL_MAJOR, 'f', k, ione, D, ione, Wloc1);
double norm_sv = LAPACKE_dlange_work( LAPACK_COL_MAJOR, 'f', k, ione, W, ione, Wloc1);
if ( myrank_mpi == 0 ) {
fprintf(stderr, "/////////////////////////////////////////////////////////////////////////\n");
}
double *ResR = (double *)malloc(k*sizeof(double));
double *ResL = (double *)malloc(k*sizeof(double));
double max_resR_err, max_resL_err;
//double *VT, *U;
//VT = (double *)malloc(n*sizeS*sizeof(double)) ;
//U = (double *)malloc(m*sizeS*sizeof(double)) ;
/* checking orthogonality U */
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &k, &k,
&alpha, &beta,
A, &i1, &i1, descA);
alpha = 1.0; beta = -1.0;
pdgemm_( "T", "N",
&k, &k, &m,
&alpha, U, &i1, &i1, descU,
U, &i1, &i1, descU,
&beta, A, &i1, &i1, descA);
double orth_U = pdlange_ ( "f", &k, &k, A, &i1, &i1, descA, Wloc1);
/* checking orthogonality VT */
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &k, &k,
&alpha, &beta,
A, &i1, &i1, descA);
alpha = 1.0; beta = -1.0;
pdgemm_( "T", "N",
&k, &k, &n,
&alpha, VT, &i1, &i1, descVT,
VT, &i1, &i1, descVT,
&beta, A, &i1, &i1, descA);
double orth_VT = pdlange_ ( "f", &k, &k, A, &i1, &i1, descA, Wloc1);
/*
* checking Berr = A-USVT
* res(ii) = norm( A*V0(:,ii)-U0(:,ii)*S0(ii,ii) );
* resL(ii) = norm(U0(:,ii)'*A-V0(:,ii)'*S0(ii,ii));
*/
/* Right Residual */
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, VT, &i1, &i1, descVT, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[0]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, U, &i1, &i1, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
ResR[0] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
int i1_i;
max_resR_err = ResR[0];
for (i = 1; i < k; i++){
i1_i = i1 + i;
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, VT, &i1, &i1_i, descVT, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[i]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, U, &i1, &i1_i, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
ResR[i] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
if ( ResR[i] > max_resR_err ) { max_resR_err = ResR[i];}
}
/* Left Residual */
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, U, &i1, &i1, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[0]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, VT, &i1, &i1, descVT, &i1, &beta, A, &i1, &i1, descA, &i1);
ResL[0] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
max_resL_err = ResL[0];
for (i = 1; i < k; i++){
i1_i = i1 + i;
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, U, &i1, &i1_i, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[i]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, VT, &i1, &i1_i, descVT, &i1, &beta, A, &i1, &i1, descA, &i1);
ResL[i] = pdlange_ ( "f", &n, &ione, A, &i1, &i1, descA, Wloc1);
if ( ResL[i] > max_resL_err ) { max_resL_err = ResL[i];}
}
if ( myrank_mpi == 0 ){
fprintf(stderr, " M \tN \tNB \tNP \tP \tQ \t#It \tCond \tTime \t\tGFlops/s \tThreshold \tNum-wanted-SV \tNum-computed-SV Acc-SV \tOrth-U \tOrth-VT ResR \tResL\n"
" %d %d \t%4d \t%4d \t%3d \t%3d \t%d %2.4e \t%7.2f \t%7.2f \t%2.4e \t%d \t\t%d \t\t%2.1e %2.1e %2.1e %2.1e %2.1e\n",
//m, n, it, cond, timer, flops/timer/1.e9, s, ktmp, sizeStmp, acc_sv/norm_sv, orth_U/m, orth_VT/n, max_resR_err/(normA*max(m,n)), max_resL_err/(normA*max(m,n)));
m, n, nb, nprocs_mpi, nprow, npcol, it, cond, elapsed_qwsvd, flops/elapsed_qwsvd/1.e9, s, ktmp, sizeStmp, acc_sv/norm_sv, orth_U/k, orth_VT/k, max_resR_err/(normA*k), max_resL_err/(normA*k));
}
}//End of checking the results
else
{
if ( myrank_mpi == 0 ){
fprintf(stderr, " M \tN \tNB \tNP \tP \tQ \t#It \tCond \tTime \t\tGFlops/s \tThreshold \tNum-wanted-SV \tNum-computed-SV Acc-SV \tOrth-U \tOrth-VT ResR \tResL\n"
" %d %d \t%4d \t%4d \t%3d \t%3d \t%d %2.4e \t%7.2f \t%7.2f \t%2.4e \t%d \t\t%d \t\t%2.1e %2.1e %2.1e %2.1e %2.1e\n",
//m, n, it, cond, timer, flops/timer/1.e9, s, ktmp, sizeStmp, acc_sv/norm_sv, orth_U/m, orth_VT/n, max_resR_err/(normA*max(m,n)), max_resL_err/(normA*max(m,n)));
m, n, nb, nprocs_mpi, nprow, npcol, it, cond, elapsed_qwsvd, flops/elapsed_qwsvd/1.e9, s, ktmp, sizeStmp, 0.0, 0.0, 0.0, 0.0, 0.0);
}
}
}//End of QDWHpartial call
/************************************************************
* End QDWHpartial
************************************************************/
/************************************************************
* Call PDGESVD
************************************************************/
if( 0 )
{
lwork_pdgesvd = -1;
Work_pdgesvd = (double *)calloc(1,sizeof(double));
pdgesvd_( "V", "V", &m, &n,
A, &i1, &i1, descA,
S,
U, &i1, &i1, descU,
VT, &i1, &i1, descVT,
Work_pdgesvd, &lwork_pdgesvd, &info);
lwork_pdgesvd = (int)Work_pdgesvd[0];
Work_pdgesvd = (double *)calloc(lwork_pdgesvd,sizeof(double));
pdlacpy_( "A", &m, &n,
Acpy, &i1, &i1, descAcpy,
A, &i1, &i1, descA );
if ( myrank_mpi == 0 ) printf("\n Start PDGESVD ... \n");
my_elapsed_pdgesvd = 0.0;
my_elapsed_pdgesvd =- MPI_Wtime();
pdgesvd_( "V", "V", &m, &n,
A, &i1, &i1, descA,
S,
U, &i1, &i1, descU,
VT, &i1, &i1, descVT,
Work_pdgesvd, &lwork_pdgesvd, &info );
my_elapsed_pdgesvd += MPI_Wtime();
MPI_Allreduce( &my_elapsed_pdgesvd, &elapsed_pdgesvd, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
sumtime_pdgesvd += elapsed_pdgesvd;
if ( elapsed_pdgesvd >= max_time_pdgesvd ) { max_time_pdgesvd = elapsed_pdgesvd;}
if ( elapsed_pdgesvd <= min_time_pdgesvd ) { min_time_pdgesvd = elapsed_pdgesvd;}
if ( myrank_mpi == 0 ) printf("\n Done PDGESVD ... \n");
long int min_Ns = min_mn * min_mn * min_mn;
long int max_Ns = max_mn * min_mn * min_mn;
if ( m!=n ){
// = (2*max_mn*min_mn^2 - 2/3*min_mn^3) + (22*min_mn^3) + 2*min_mn^3
flops_pdgesvd += 2.*(double)max_Ns -
2./3.*(double)min_Ns +
22.*(double)min_Ns +
2.*(double)min_Ns;
}
else{
flops_pdgesvd += 22.*(double)min_Ns;
}
/*
* Check the results
*/
if ( check )
{
if ( myrank_mpi == 0 ) {
fprintf(stderr, "\n PDGESVD \n");
fprintf(stderr, "/////////////////////////////////////////////////////////////////////////\n");
}
dlasrt_( "D", &n, W, &info );
dlasrt_( "D", &n, S, &info );
for (i = 0; i < k; i++){
D[i] = W[i] - S[i];
}
alpha = 1.0; beta = -1.0;
int ione = 1;
double acc_sv = LAPACKE_dlange_work( LAPACK_COL_MAJOR, 'f', k, ione, D, ione, Wloc1);
double norm_sv = LAPACKE_dlange_work( LAPACK_COL_MAJOR, 'f', k, ione, W, ione, Wloc1);
if ( myrank_mpi == 0 ) {
fprintf(stderr, "/////////////////////////////////////////////////////////////////////////\n");
}
double *ResR = (double *)malloc(k*sizeof(double));
double *ResL = (double *)malloc(k*sizeof(double));
double max_resR_err, max_resL_err;
//double *VT, *U;
//VT = (double *)malloc(n*sizeS*sizeof(double)) ;
//U = (double *)malloc(m*sizeS*sizeof(double)) ;
/* checking orthogonality U */
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &k, &k,
&alpha, &beta,
A, &i1, &i1, descA);
alpha = 1.0; beta = -1.0;
pdgemm_( "T", "N",
&k, &k, &m,
&alpha, U, &i1, &i1, descU,
U, &i1, &i1, descU,
&beta, A, &i1, &i1, descA);
double orth_U = pdlange_ ( "f", &k, &k, A, &i1, &i1, descA, Wloc1);
/* checking orthogonality VT */
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &k, &k,
&alpha, &beta,
A, &i1, &i1, descA);
alpha = 1.0; beta = -1.0;
pdgemm_( "T", "N",
&k, &k, &n,
&alpha, VT, &i1, &i1, descVT,
VT, &i1, &i1, descVT,
&beta, A, &i1, &i1, descA);
double orth_VT = pdlange_ ( "f", &k, &k, A, &i1, &i1, descA, Wloc1);
/*
* checking Berr = A-USVT
* res(ii) = norm( A*V0(:,ii)-U0(:,ii)*S0(ii,ii) );
* resL(ii) = norm(U0(:,ii)'*A-V0(:,ii)'*S0(ii,ii));
*/
/* Right Residual */
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, VT, &i1, &i1, descVT, &n, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[0]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, U, &i1, &i1, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
ResR[0] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
int i1_i;
max_resR_err = ResR[0];
for (i = 1; i < k; i++){
i1_i = i1 + i;
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, VT, &i1_i, &i1, descVT, &n, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[i]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, U, &i1, &i1_i, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
ResR[i] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
if ( ResR[i] > max_resR_err ) { max_resR_err = ResR[i];}
}
/* Left Residual */
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, U, &i1, &i1, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[0]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, VT, &i1, &i1, descVT, &n, &beta, A, &i1, &i1, descA, &i1);
ResL[0] = pdlange_ ( "f", &m, &ione, A, &i1, &i1, descA, Wloc1);
max_resL_err = ResL[0];
for (i = 1; i < k; i++){
i1_i = i1 + i;
alpha = 1.0; beta = 0.0;
pdgemv_ ("N", &m, &n, &alpha, Acpy, &i1, &i1, descAcpy, U, &i1, &i1_i, descU, &i1, &beta, A, &i1, &i1, descA, &i1);
alpha = 0.0; beta = 1.0;
pdlaset_( "G", &m, &n,
&alpha, &beta,
pinv, &i1, &i1, descpinv);
alpha = -S[i]; beta = 1.0;
pdgemv_ ("N", &m, &n, &alpha, pinv, &i1, &i1, descpinv, VT, &i1_i, &i1, descVT, &n, &beta, A, &i1, &i1, descA, &i1);
ResL[i] = pdlange_ ( "f", &n, &ione, A, &i1, &i1, descA, Wloc1);
if ( ResL[i] > max_resL_err ) { max_resL_err = ResL[i];}
}
if ( myrank_mpi == 0 ){
fprintf(stderr, " M \tN \tNB \tNP \tP \tQ \t#It \tCond \tTime \t\tGFlops/s \tThreshold \tNum-wanted-SV \tNum-computed-SV Acc-SV \tOrth-U \tOrth-VT ResR \tResL\n"
" %d %d \t%4d \t%4d \t%3d \t%3d \t%d %2.4e \t%7.2f \t%7.2f \t%2.4e \t%d \t\t%d \t\t%2.1e %2.1e %2.1e %2.1e %2.1e\n",
m, n, nb, nprocs_mpi, nprow, npcol, 0, cond, elapsed_pdgesvd, flops_pdgesvd/elapsed_pdgesvd/1.e9, 0.0, ktmp, sizeStmp, acc_sv/norm_sv, orth_U/k, orth_VT/k, max_resR_err/(normA*k), max_resL_err/(normA*k));
}
}//End of checking the results
else
{
if ( myrank_mpi == 0 ){
fprintf(stderr, " M \tN \tNB \tNP \tP \tQ \t#It \tCond \tTime \t\tGFlops/s \tThreshold \tNum-wanted-SV \tNum-computed-SV Acc-SV \tOrth-U \tOrth-VT ResR \tResL\n"
" %d %d \t%4d \t%4d \t%3d \t%3d \t%d %2.4e \t%7.2f \t%7.2f \t%2.4e \t%d \t\t%d \t\t%2.1e %2.1e %2.1e %2.1e %2.1e\n",
m, n, nb, nprocs_mpi, nprow, npcol, 0, cond, elapsed_pdgesvd, flops_pdgesvd/elapsed_pdgesvd/1.e9, 0.0, ktmp, sizeStmp, 0.0, 0.0, 0.0, 0.0, 0.0);
}
}
}//End of PDGESVD call
/************************************************************
* End PDGESVD
************************************************************/
} // loop over iteration
free(A); free(B); free(U); free(VT); free(pinv); free(S); free(Acpy); free(D); free(W);
} // End loop over range
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Range loop ends\n");
blacs_gridexit_( &i0 );
MPI_Finalize();
if ( verbose & myrank_mpi == 0 ) fprintf(stderr, "Program ends...\n");
return 0;
}
int
get_range(char *range, int *start_p, int *stop_p, int *step_p) {
char *s, *s1, buf[21];
int colon_count, copy_len, nbuf=20, n;
int start=1000, stop=10000, step=1000;
colon_count = 0;
for (s = strchr( range, ':'); s; s = strchr( s+1, ':'))
colon_count++;
if (colon_count == 0) { /* No colon in range. */
if (sscanf( range, "%d", &start ) < 1 || start < 1)
return -1;
step = start / 10;
if (step < 1) step = 1;
stop = start + 10 * step;
} else if (colon_count == 1) { /* One colon in range.*/
/* First, get the second number (after colon): the stop value. */
s = strchr( range, ':' );
if (sscanf( s+1, "%d", &stop ) < 1 || stop < 1)
return -1;
/* Next, get the first number (before colon): the start value. */
n = s - range;
copy_len = n > nbuf ? nbuf : n;
strncpy( buf, range, copy_len );
buf[copy_len] = 0;
if (sscanf( buf, "%d", &start ) < 1 || start > stop || start < 1)
return -1;
/* Let's have 10 steps or less. */
step = (stop - start) / 10;
if (step < 1)
step = 1;
} else if (colon_count == 2) { /* Two colons in range. */
/* First, get the first number (before the first colon): the start value. */
s = strchr( range, ':' );
n = s - range;
copy_len = n > nbuf ? nbuf : n;
strncpy( buf, range, copy_len );
buf[copy_len] = 0;
if (sscanf( buf, "%d", &start ) < 1 || start < 1)
return -1;
/* Next, get the second number (after the first colon): the stop value. */
s1 = strchr( s+1, ':' );
n = s1 - (s + 1);
copy_len = n > nbuf ? nbuf : n;
strncpy( buf, s+1, copy_len );
buf[copy_len] = 0;
if (sscanf( buf, "%d", &stop ) < 1 || stop < start)
return -1;
/* Finally, get the third number (after the second colon): the step value. */
if (sscanf( s1+1, "%d", &step ) < 1 || step < 1)
return -1;
} else
return -1;
*start_p = start;
*stop_p = stop;
*step_p = step;
return 0;
}