-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpsymbfact.c
5233 lines (4825 loc) · 162 KB
/
psymbfact.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
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*! \file
Copyright (c) 2003, The Regents of the University of California, through
Lawrence Berkeley National Laboratory (subject to receipt of any required
approvals from U.S. Dept. of Energy)
All rights reserved.
The source code is distributed under BSD license, see the file License.txt
at the top-level directory.
*/
/*! @file
* \brief Implements parallel symbolic factorization
*
* <pre>
* -- Parallel symbolic factorization routine (version 2.3) --
* Lawrence Berkeley National Lab, Univ. of California Berkeley - July 2003
* INRIA France - January 2004
* Laura Grigori
*
* November 1, 2007
* Feburary 20, 2008
* October 15, 2008
* January 28, 2018
*
* The function symbfact_dist implements the parallel symbolic factorization
* algorithm described in the paper:
*
* Parallel Symbolic Factorization for Sparse LU with Static Pivoting,
* Laura Grigori, James W. Demmel and Xiaoye S. Li,
* Pages 1289-1314, SIAM Journal on Scientific Computing, Volume 29, Issue 3.
* </pre>
*/
/* limits.h: the largest positive integer (INT_MAX) (LONG_MAX) */
#include <limits.h>
#include <math.h>
#include "superlu_ddefs.h"
#include "psymbfact.h"
/*
* Internal protypes
*/
static int_t *
intMalloc_symbfact(int_t );
static int_t *
intCalloc_symbfact(int_t );
static int_t
initParmsAndStats
(psymbfact_stat_t *PS);
static void
estimate_memUsage
(int_t, int, superlu_dist_mem_usage_t *, float *, float *,
Pslu_freeable_t *, Llu_symbfact_t *,
vtcsInfo_symbfact_t *, comm_symbfact_t *, psymbfact_stat_t *);
static void
symbfact_free
(int, int, Llu_symbfact_t *, vtcsInfo_symbfact_t *, comm_symbfact_t *);
static int_t
denseSep_symbfact
(int , int_t, int, int, int, int_t *, int_t *, int,
int, int, int_t, int_t, int_t *, int_t *, int_t *,
int_t *, int_t *, MPI_Comm, MPI_Comm *, Llu_symbfact_t *,
Pslu_freeable_t *_freeable, vtcsInfo_symbfact_t *,
comm_symbfact_t *, psymbfact_stat_t * );
static int_t
dnsUpSeps_symbfact
(int_t, int, int, int, int, int_t *, int_t *, int_t,
Llu_symbfact_t *, Pslu_freeable_t *, vtcsInfo_symbfact_t *,
comm_symbfact_t *, psymbfact_stat_t *, int_t *, int_t *, int_t *);
static void
intraLvl_symbfact
(SuperMatrix *, int, int, int, int, int, int_t *, int_t *, int,
int, int_t, int_t, Pslu_freeable_t *, Llu_symbfact_t *, vtcsInfo_symbfact_t *,
comm_symbfact_t *, psymbfact_stat_t *, int_t *, int_t *, int_t *, int_t *,
int_t *, int_t *, int_t *, MPI_Comm, MPI_Comm *);
static void
initLvl_symbfact
(int_t, int, int_t, int_t, Pslu_freeable_t *,
Llu_symbfact_t *, vtcsInfo_symbfact_t *, psymbfact_stat_t *, MPI_Comm,
int_t *, int_t, int_t);
static void
createComm (int, int, MPI_Comm *, MPI_Comm *);
static void
freeComm (int, int, MPI_Comm *, MPI_Comm *);
static void
domain_symbfact
(SuperMatrix *, int, int, int, int, int, int_t *, int_t *,
int_t, int_t, Pslu_freeable_t *, Llu_symbfact_t *, vtcsInfo_symbfact_t *,
comm_symbfact_t *, psymbfact_stat_t *, int_t *, int_t *, int_t *, int_t *,
int_t *, int_t *, int_t *);
static float
allocPrune_domain
(int_t, int_t, Llu_symbfact_t *, vtcsInfo_symbfact_t *, psymbfact_stat_t *);
static float
allocPrune_lvl
(Llu_symbfact_t *, vtcsInfo_symbfact_t *, psymbfact_stat_t *);
static int
symbfact_alloc
(int_t, int, Pslu_freeable_t *, Llu_symbfact_t *,
vtcsInfo_symbfact_t *, comm_symbfact_t *, psymbfact_stat_t *);
static float
symbfact_mapVtcs
(int, int, int, SuperMatrix *, int_t *, int_t *,
Pslu_freeable_t *, vtcsInfo_symbfact_t *, int_t *, int_t, psymbfact_stat_t *);
static void
symbfact_distributeMatrix
(int, int, int, SuperMatrix *, int_t *, int_t *, matrix_symbfact_t *,
Pslu_freeable_t *, vtcsInfo_symbfact_t *, int_t *, MPI_Comm *);
static int_t
interLvl_symbfact
(SuperMatrix *, int, int, int, int, int, int, int,
int_t *, int_t *, int_t *, int_t *, int_t *, int_t *, int_t *,
Llu_symbfact_t *, Pslu_freeable_t*, comm_symbfact_t *, vtcsInfo_symbfact_t *,
psymbfact_stat_t *, MPI_Comm, MPI_Comm *);
static float
cntsVtcs
(int_t, int, int, Pslu_freeable_t *, Llu_symbfact_t *, vtcsInfo_symbfact_t *,
int_t *, int_t *, int_t *, psymbfact_stat_t *, MPI_Comm *);
/************************************************************************/
float symbfact_dist
/************************************************************************/
(
int nprocs_num, /* Input - no of processors */
int nprocs_symb, /* Input - no of processors for the symbolic
factorization */
SuperMatrix *A, /* Input - distributed input matrix */
int_t *perm_c, /* Input - column permutation */
int_t *perm_r, /* Input - row permutation */
int_t *sizes, /* Input - sizes of each node in the separator tree */
int_t *fstVtxSep, /* Input - first vertex of each node in the tree */
Pslu_freeable_t *Pslu_freeable, /* Output - local L and U structure,
global to local indexing information */
MPI_Comm *num_comm, /* Input - communicator for numerical factorization */
MPI_Comm *symb_comm, /* Input - communicator for symbolic factorization */
superlu_dist_mem_usage_t *symb_mem_usage
)
{
/*! \brief
*
* <pre>
* Purpose
* =======
* symbfact_dist() performs symbolic factorization of matrix A suitable
* for performing the supernodal Gaussian elimination with no pivoting (GEPP).
* This routine computes the structure of one column of L and one row of U
* at a time. It uses:
* o distributed input matrix
* o supernodes
* o symmetric structure pruning
*
*
* Arguments
* =========
*
* nprocs_num (input) int
* Number of processors SuperLU_DIST is executed on, and the input
* matrix is distributed on.
*
* nprocs_symb (input) int
* Number of processors on which the symbolic factorization is
* performed. It is equal to the number of independent domains
* idenfied in the graph partitioning algorithm executed
* previously and has to be a power of 2. It corresponds to
* number of leaves in the separator tree.
*
* A (input) SuperMatrix*
* Matrix A in A*X=B, of dimension (A->nrow, A->ncol). The
* number of the linear equations is A->nrow. Matrix A is
* distributed in NRformat_loc format.
* Matrix A is not yet permuted by perm_c.
*
* perm_c (input) int_t*
* Column permutation vector of size A->ncol, which defines the
* permutation matrix Pc; perm_c[i] = j means column i of A is
* in position j in A*Pc.
*
* perm_r (input) int_t*
* Row permutation vector of size A->nrow, which defines the
* permutation matrix Pr; perm_r[i] = j means column i of A is
* in position j in Pr*A.
*
* sizes (input) int_t*
* Contains the number of vertices in each separator.
*
* fstVtxSep (input) int_t*
* Contains first vertex for each separator.
*
* Pslu_freeable (output) Pslu_freeable_t*
* Returns the local L and U structure, and global to local
* information on the indexing of the vertices. Contains all
* the information necessary for performing the data
* distribution towards the numeric factorization.
*
* num_comm (input) MPI_Comm*
* Communicator for numerical factorization
*
* symb_comm (input) MPI_Comm*
* Communicator for symbolic factorization
*
* symb_mem_usage (input) superlu_dist_mem_usage_t *
* Statistics on memory usage.
*
* Return value
* ============
* < 0, number of bytes allocated on return from the symbolic factorization.
* > 0, number of bytes allocated when out of memory.
*
* Sketch of the algorithm
* =======================
*
* Distrbute the vertices on the processors using a subtree to
* subcube algorithm.
*
* Redistribute the structure of the input matrix A according to the
* subtree to subcube computed previously for the symbolic
* factorization routine. This implies in particular a distribution
* from nprocs_num processors to nprocs_symb processors.
*
* Perform symbolic factorization guided by the separator tree provided by
* a graph partitioning algorithm. The symbolic factorization uses a
* combined left-looking, right-looking approach.
* </pre>
*/
NRformat_loc *Astore;
int iam, szSep, fstP, lstP, npNode, nlvls, lvl, p, iSep, jSep;
int iinfo; /* return code */
int_t m, n;
int_t nextl, nextu, neltsZr, neltsTotal, nsuper_loc, szLGr, szUGr;
int_t ind_blk, nsuper, vtx, min_mn, szsn;
long long int nnzL, nnzU, nnzLU;
float stat_loc[23], stat_glob[23], mem_glob[15];
Llu_symbfact_t Llu_symbfact; /* local L and U and pruned L and U data structures */
vtcsInfo_symbfact_t VInfo; /* local information on number of blocks,
number of vertices in a block etc */
matrix_symbfact_t AS; /* temporary storage for the input matrix after redistribution */
comm_symbfact_t CS; /* information on communication */
/* relaxation parameters (for future release) and
statistics collected during the symbolic factorization */
psymbfact_stat_t PS;
/* temp array of size n, used as a marker by the subroutines */
int_t *tempArray;
int_t i, j, k;
int_t fstVtx, lstVtx, mark, fstVtx_lid, vtx_lid, maxNvtcsPProc;
int_t nnz_asup_loc, nnz_ainf_loc, fill_rcmd;
float totalMemLU, overestimMem;
MPI_Comm *commLvls;
/* maximum block size */
int_t maxSzBlk;
float flinfo;
#if ( PRNTlevel >= 1)
float stat_msgs_l[10], stat_msgs_g[10];
#endif
#if ( PROFlevel>=1 )
double t, t_symbFact[3], t_symbFact_loc[3];
double *time_lvlsT, *time_lvls, t1, t2, time_lvlsg[9];
#endif
/* Initialization */
MPI_Comm_rank ((*num_comm), &iam);
commLvls = NULL;
#if ( DEBUGlevel>=1 )
CHECK_MALLOC(iam, "Enter psymbfact()");
#endif
initParmsAndStats (&PS);
if (nprocs_symb != 1) {
if (!(commLvls = (MPI_Comm *) SUPERLU_MALLOC(2*nprocs_symb*sizeof(MPI_Comm)))) {
fprintf (stderr, "Malloc fails for commLvls[].");
return (PS.allocMem);
}
PS.allocMem += 2 * nprocs_symb * sizeof(MPI_Comm);
}
nlvls = (int) LOG2( nprocs_num ) + 1;
#if ( PROFlevel>=1 )
time_lvlsT = (double *) SUPERLU_MALLOC(3*nprocs_symb*(nlvls+1)
* sizeof(double));
time_lvls = (double *) SUPERLU_MALLOC(3*(nlvls+1) * sizeof(double));
if (!time_lvls || !time_lvlsT) {
fprintf (stderr, "Malloc fails for time_lvls[].");
return (PS.allocMem);
}
PS.allocMem += (3*nprocs_symb*(nlvls+1) + 3*(nlvls+1)) * sizeof(double);
#endif
VInfo.xlsub_nextLvl = 0;
VInfo.xusub_nextLvl = 0;
VInfo.maxSzBlk = sp_ienv_dist(3);
maxSzBlk = VInfo.maxSzBlk;
mark = EMPTY;
nsuper_loc = 0;
nextl = 0; nextu = 0;
neltsZr = 0; neltsTotal = 0;
m = A->nrow;
n = A->ncol;
min_mn = SUPERLU_MIN( m, n );
if (!(tempArray = intMalloc_symbfact(n))) {
fprintf (stderr, "Malloc fails for tempArray[].\n");
return (PS.allocMem);
}
PS.allocMem += n * sizeof(int_t);
#if ( PROFlevel>=1 )
t = SuperLU_timer_();
#endif
/* Distribute vertices on processors */
if ((flinfo =
symbfact_mapVtcs (iam, nprocs_num, nprocs_symb, A, fstVtxSep, sizes,
Pslu_freeable, &VInfo, tempArray, maxSzBlk, &PS)) > 0)
return (flinfo);
maxNvtcsPProc = Pslu_freeable->maxNvtcsPProc;
/* Redistribute matrix A on processors following the distribution found
in symbfact_mapVtcs. Store the redistributed A temporarily into AS */
symbfact_distributeMatrix (iam, nprocs_num, nprocs_symb, A,
perm_c, perm_r, &AS,
Pslu_freeable, &VInfo, tempArray, num_comm);
/* THE REST OF THE SYMBOLIC FACTORIZATION IS EXECUTED ONLY BY NPROCS_SYMB
PROCESSORS */
if ( iam < nprocs_symb ) {
#if ( PROFlevel>=1 )
t_symbFact_loc[0] = SuperLU_timer_() - t;
t = SuperLU_timer_();
t_symbFact_loc[1] = t;
#endif
/* Allocate storage common to the symbolic factor routines */
if (iinfo = symbfact_alloc (n, nprocs_symb, Pslu_freeable,
&Llu_symbfact, &VInfo, &CS, &PS))
return (PS.allocMem);
/* Copy the redistributed input matrix AS at the end of the memory buffer
allocated to store L and U. That is, copy (AS.x_ainf, AS.ind_ainf) in
(xlsub, lsub), (AS.x_asup, AS.ind_asup) in (xusub, usub). Free the
memory used to store the input matrix */
nnz_ainf_loc = VInfo.nnz_ainf_loc;
nnz_asup_loc = VInfo.nnz_asup_loc;
j = Llu_symbfact.szUsub - VInfo.nnz_asup_loc;
k = Llu_symbfact.szLsub - VInfo.nnz_ainf_loc;
for (i = 0; i <= VInfo.nvtcs_loc; i++) {
Llu_symbfact.xusub[i] = AS.x_asup[i] + j;
Llu_symbfact.xlsub[i] = AS.x_ainf[i] + k;
}
for (i = 0; i < VInfo.nnz_asup_loc; i++, j++)
Llu_symbfact.usub[j] = AS.ind_asup[i];
for (i = 0; i < VInfo.nnz_ainf_loc; i++, k++)
Llu_symbfact.lsub[k] = AS.ind_ainf[i];
SUPERLU_FREE( AS.x_ainf );
SUPERLU_FREE( AS.x_asup );
SUPERLU_FREE( AS.ind_ainf );
SUPERLU_FREE( AS.ind_asup );
if (nprocs_symb != 1) {
createComm (iam, nprocs_symb, commLvls, symb_comm);
#if ( PROFlevel>=1 )
t_symbFact_loc[2] = SuperLU_timer_();
#endif
if ((flinfo = cntsVtcs (n, iam, nprocs_symb, Pslu_freeable, &Llu_symbfact,
&VInfo, tempArray, fstVtxSep, sizes, &PS, commLvls)) > 0)
return (flinfo);
#if ( PROFlevel>=1 )
t_symbFact_loc[2] = SuperLU_timer_() - t_symbFact_loc[2];
#endif
}
/* set to EMPTY marker[] array */
for (i = 0; i < n; i++)
tempArray[i] = EMPTY;
szSep = nprocs_symb;
iSep = 0;
lvl = 0;
while (szSep >= 1) {
/* for each level in the separator tree */
npNode = nprocs_symb / szSep;
fstP = 0;
/* for each node in the level */
for (jSep = iSep; jSep < iSep + szSep; jSep++) {
fstVtx = fstVtxSep[jSep];
lstVtx = fstVtx + sizes[jSep];
/* if this is the first level */
if (szSep == nprocs_symb) {
/* compute symbolic factorization for my domain */
if (fstP == iam) {
/* allocate storage for the pruned structures */
#if ( PROFlevel>=1 )
t1 = SuperLU_timer_();
#endif
if ((flinfo = allocPrune_domain (fstVtx, lstVtx,
&Llu_symbfact, &VInfo, &PS)) > 0)
return (flinfo);
if (fstVtx < lstVtx)
VInfo.fstVtx_nextLvl = VInfo.begEndBlks_loc[2];
domain_symbfact
(A, iam, lvl, szSep, iSep, jSep, sizes, fstVtxSep, fstVtx, lstVtx,
Pslu_freeable, &Llu_symbfact, &VInfo, &CS, &PS, tempArray,
&mark, &nextl, &nextu, &neltsZr, &neltsTotal, &nsuper_loc);
PS.estimLSz = nextl;
PS.estimUSz = nextu;
if (nprocs_symb != 1)
if((flinfo = allocPrune_lvl (&Llu_symbfact, &VInfo, &PS)) > 0)
return (flinfo);
#if ( PROFlevel>=1 )
t2 = SuperLU_timer_();
time_lvls[lvl] = 0.; time_lvls[lvl+1] = 0.;
time_lvls[lvl + 2] = t2 - t1;
#endif
}
}
else {
lstP = fstP + npNode;
if (fstP <= iam && iam < lstP) {
#if ( PROFlevel>=1 )
t1 = SuperLU_timer_();
#endif
if (VInfo.filledSep != FILLED_SEPS)
initLvl_symbfact(n, iam, fstVtx, lstVtx,
Pslu_freeable, &Llu_symbfact, &VInfo, &PS, commLvls[jSep],
tempArray, nextl, nextu);
#if ( PROFlevel>=1 )
t2 = SuperLU_timer_();
time_lvls[3*lvl] = t2 - t1;
#endif
interLvl_symbfact (A, iam, lvl, szSep, fstP, lstP,
iSep, jSep, sizes, fstVtxSep,
&nextl, &nextu, &nsuper_loc, &mark, tempArray,
&Llu_symbfact, Pslu_freeable, &CS, &VInfo, &PS,
commLvls[jSep], symb_comm);
#if ( PROFlevel>=1 )
t1 = SuperLU_timer_();
time_lvls[3*lvl+1] = t1 - t2;
#endif
if (VInfo.filledSep != FILLED_SEPS)
intraLvl_symbfact
(A, iam, lvl, szSep, iSep, jSep, sizes, fstVtxSep, fstP, lstP,
fstVtx, lstVtx, Pslu_freeable, &Llu_symbfact, &VInfo, &CS, &PS,
tempArray, &mark, &nextl, &nextu, &neltsZr, &neltsTotal,
&nsuper_loc, commLvls[jSep], symb_comm);
#if ( PROFlevel>=1 )
t2 = SuperLU_timer_();
time_lvls[3*lvl+2] = t2 - t1;
#endif
}
}
fstP += npNode;
}
iSep += szSep;
szSep = szSep / 2;
lvl ++;
}
SUPERLU_FREE( tempArray );
/* Set up global information and collect statistics */
if (PS.maxSzLPr < Llu_symbfact.indLsubPr)
PS.maxSzLPr = Llu_symbfact.indLsubPr;
if (PS.maxSzUPr < Llu_symbfact.indUsubPr)
PS.maxSzUPr = Llu_symbfact.indUsubPr;
Llu_symbfact.xlsub[VInfo.nvtcs_loc] = nextl;
Llu_symbfact.xusub[VInfo.nvtcs_loc] = nextu;
fill_rcmd = SUPERLU_MAX( nextl / (nnz_ainf_loc+1), nextu / (nnz_asup_loc+1)) + 1;
Pslu_freeable->xsup_beg_loc = intMalloc_dist (nsuper_loc+1);
Pslu_freeable->xsup_end_loc = intMalloc_dist (nsuper_loc+1);
if (!Pslu_freeable->xsup_beg_loc || !Pslu_freeable->xsup_end_loc) {
fprintf (stderr, "Malloc fails for xsup_beg_loc, xsup_end_loc.");
return (PS.allocMem);
}
PS.allocMem += 2 * (nsuper_loc+1) * sizeof(int_t);
maxNvtcsPProc = Pslu_freeable->maxNvtcsPProc;
nnzL = 0; nnzU = 0;
i = 0;
nsuper = 0;
ind_blk = 0;
for (ind_blk = 0; ind_blk < VInfo.nblks_loc; ind_blk ++) {
fstVtx = VInfo.begEndBlks_loc[2 * ind_blk];
lstVtx = VInfo.begEndBlks_loc[2 * ind_blk + 1];
fstVtx_lid = LOCAL_IND( Pslu_freeable->globToLoc[fstVtx] );
nsuper = Pslu_freeable->supno_loc[fstVtx_lid];
Pslu_freeable->xsup_beg_loc[nsuper] = fstVtx;
szsn = 1;
if (LONG_MAX - nnzL <= Llu_symbfact.xlsub[fstVtx_lid + 1] -
Llu_symbfact.xlsub[fstVtx_lid])
printf ("PE[%d] ERR nnzL %lld\n", iam, nnzL);
if (LONG_MAX - nnzU <= Llu_symbfact.xusub[fstVtx_lid + 1] -
Llu_symbfact.xusub[fstVtx_lid])
printf ("PE[%d] ERR nnzU %lld\n", iam, nnzU);
j = Llu_symbfact.xlsub[fstVtx_lid + 1] - Llu_symbfact.xlsub[fstVtx_lid];
k = Llu_symbfact.xusub[fstVtx_lid + 1] - Llu_symbfact.xusub[fstVtx_lid];
nnzL += j;
nnzU += k;
for (vtx = fstVtx + 1, vtx_lid = fstVtx_lid + 1;
vtx < lstVtx; vtx++, vtx_lid ++) {
if (Pslu_freeable->supno_loc[vtx_lid] != nsuper) {
nsuper = Pslu_freeable->supno_loc[vtx_lid];
Pslu_freeable->xsup_end_loc[nsuper-1] = vtx;
Pslu_freeable->xsup_beg_loc[nsuper] = vtx;
szsn = 1;
j = Llu_symbfact.xlsub[vtx_lid + 1] - Llu_symbfact.xlsub[vtx_lid];
k = Llu_symbfact.xusub[vtx_lid + 1] - Llu_symbfact.xusub[vtx_lid];
}
else {
szsn ++;
}
nnzL += j - szsn + 1;
nnzU += k - szsn + 1;
}
Pslu_freeable->xsup_end_loc[nsuper] = lstVtx;
}
Pslu_freeable->supno_loc[VInfo.nvtcs_loc] = nsuper_loc;
Pslu_freeable->nvtcs_loc = VInfo.nvtcs_loc;
/* set up xsup data */
Pslu_freeable->lsub = Llu_symbfact.lsub;
Pslu_freeable->xlsub = Llu_symbfact.xlsub;
Pslu_freeable->usub = Llu_symbfact.usub;
Pslu_freeable->xusub = Llu_symbfact.xusub;
Pslu_freeable->szLsub = Llu_symbfact.szLsub;
Pslu_freeable->szUsub = Llu_symbfact.szUsub;
#if ( PROFlevel>=1 )
t_symbFact_loc[1] = SuperLU_timer_() - t_symbFact_loc[1];
#endif
#if ( PRNTlevel>=1 )
estimate_memUsage (n, iam, symb_mem_usage,
&totalMemLU, &overestimMem,
Pslu_freeable, &Llu_symbfact, &VInfo, &CS, &PS);
stat_loc[0] = (float) nnzL;
stat_loc[1] = (float) nnzU;
stat_loc[2] = (float) nsuper_loc;
stat_loc[3] = (float) Pslu_freeable->xlsub[VInfo.nvtcs_loc];
stat_loc[4] = (float) Pslu_freeable->xusub[VInfo.nvtcs_loc];
stat_loc[5] = totalMemLU;
stat_loc[6] = overestimMem;
stat_loc[7] = totalMemLU - overestimMem;
stat_loc[8] = (float) PS.maxSzBuf;
stat_loc[9] = (float) PS.nDnsUpSeps;
stat_loc[10] = (float) PS.nDnsCurSep;
stat_loc[11] = (float) (Llu_symbfact.no_expand + Llu_symbfact.no_expcp +
Llu_symbfact.no_expand_pr);
stat_loc[12] = (float) Llu_symbfact.no_expand;
stat_loc[13] = (float) Llu_symbfact.no_expcp;
stat_loc[14] = (float) Llu_symbfact.no_expand_pr;
stat_loc[15] = (float) fill_rcmd;
stat_loc[16] = PS.nops;
stat_loc[17] = PS.fill_pelt[1];
stat_loc[18] = PS.fill_pelt[4];
stat_loc[19] = PS.fill_pelt[0];
stat_loc[20] = PS.fill_pelt[2];
stat_loc[21] = PS.fill_pelt[3];
stat_loc[22] = PS.fill_pelt[5];
MPI_Allreduce (stat_loc, stat_glob, 23, MPI_FLOAT,
MPI_SUM, (*symb_comm));
MPI_Allreduce (&(stat_loc[5]), mem_glob, 14, MPI_FLOAT,
MPI_MAX, (*symb_comm));
fill_rcmd = (int_t) mem_glob[10];
PS.fill_pelt[0] = stat_glob[19];
PS.fill_pelt[1] = mem_glob[12];
PS.fill_pelt[2] = stat_glob[20];
PS.fill_pelt[3] = stat_glob[21];
PS.fill_pelt[4] = mem_glob[13];
PS.fill_pelt[5] = stat_glob[22];
if (PS.fill_pelt[2] == 0.) PS.fill_pelt[2] = 1.;
if (PS.fill_pelt[5] == 0.) PS.fill_pelt[5] = 1.;
#if ( PROFlevel>=1 )
MPI_Reduce (t_symbFact_loc, t_symbFact, 3, MPI_DOUBLE,
MPI_MAX, 0, (*symb_comm));
MPI_Gather (time_lvls, 3 * nlvls, MPI_DOUBLE,
time_lvlsT, 3 * nlvls , MPI_DOUBLE,
0, (*symb_comm));
#endif
stat_msgs_l[0] = (float) PS.maxsz_msgSnd;
stat_msgs_l[1] = (float) PS.maxsz_msgSnd;
if (PS.maxsz_msgSnd < PS.maxsz_msgCol)
stat_msgs_l[1] = PS.maxsz_msgCol;
stat_msgs_l[2] = PS.no_shmSnd + PS.no_msgsSnd +
PS.no_shmRcvd + PS.no_msgsRcvd;
stat_msgs_l[3] = stat_msgs_l[2] + PS.no_msgsCol;
stat_msgs_l[4] = stat_msgs_l[2];
stat_msgs_l[5] = stat_msgs_l[3];
stat_msgs_l[6] = PS.no_msgsSnd;
stat_msgs_l[7] = PS.no_msgsSnd + PS.no_msgsCol;
stat_msgs_l[8] = PS.sz_msgsSnd;
stat_msgs_l[9] = PS.sz_msgsSnd + PS.sz_msgsCol;
MPI_Reduce (stat_msgs_l, stat_msgs_g, 4, MPI_FLOAT,
MPI_MAX, 0, (*symb_comm));
MPI_Reduce (&(stat_msgs_l[4]), &(stat_msgs_g[4]), 6, MPI_FLOAT,
MPI_SUM, 0, (*symb_comm));
if (stat_msgs_g[6] == 0) stat_msgs_g[6] = 1;
if (stat_msgs_g[7] == 0) stat_msgs_g[7] = 1;
Pslu_freeable->nnzLU=(long long) stat_glob[0]+(long long) stat_glob[1];
if (!iam) {
nnzL = (long long) stat_glob[0]; nnzU = (long long) stat_glob[1];
nsuper = (int_t) stat_glob[2];
szLGr = (int_t) stat_glob[3]; szUGr = (int_t) stat_glob[4];
printf("\tMax szBlk %ld\n", (long) VInfo.maxSzBlk);
#if ( PRNTlevel>=2 )
printf("\t relax_gen %.2f, relax_curSep %.2f, relax_seps %.2f\n",
PS.relax_gen, PS.relax_curSep, PS.relax_seps);
#endif
printf("LONG_MAX %ld\n", LONG_MAX);
printf("\tParameters: fill mem %ld fill pelt %ld\n",
(long) sp_ienv_dist(6), (long) PS.fill_par);
printf("\tNonzeros in L %lld\n", nnzL);
printf("\tNonzeros in U %lld\n", nnzU);
nnzLU = nnzL + nnzU;
printf("\tnonzeros in L+U-I %lld\n", nnzLU);
printf("\tNo of supers %ld\n", (long) nsuper);
printf("\tSize of G(L) %ld\n", (long) szLGr);
printf("\tSize of G(U) %ld\n", (long) szUGr);
printf("\tSize of G(L+U) %ld\n", (long) szLGr+szUGr);
printf("\tParSYMBfact (MB) :\tL\\U MAX %.2f\tAVG %.2f\n",
mem_glob[0]*1e-6,
stat_glob[5]/nprocs_symb*1e-6);
#if ( PRNTlevel>=2 )
printf("\tRL overestim (MB):\tL\\U MAX %.2f\tAVG %.2f\n",
mem_glob[1]*1e-6,
stat_glob[6]/nprocs_symb*1e-6);
printf("\tsnd/rcv buffers (MB):\tL\\U MAX %.2f\tAVG %.2f\n",
mem_glob[3]*1e-6,
stat_glob[8]/nprocs_symb*1e-6);
printf("\tSYMBfact 2*n+4*nvtcs_loc+2*maxNvtcsNds_loc:\tL\\U %.2f\n",
(float) (2 * n * sizeof(int_t)) *1e-6);
printf("\tint_t %d, int %d, long int %d, short %d, float %d, double %d\n",
sizeof(int_t), sizeof(int), sizeof(long int), sizeof(short), sizeof(float),
sizeof(double));
printf("\tDNS ALLSEPS:\t MAX %d\tAVG %.2f\n",
(int_t) mem_glob[4], stat_glob[9]/nprocs_symb);
printf("\tDNS CURSEP:\t MAX %d\tAVG %.2f\n\n",
(int_t) mem_glob[5], stat_glob[10]/nprocs_symb);
printf("\t MAX FILL Mem(L+U) / Mem(A) per processor %ld\n", fill_rcmd);
printf("\t Per elt MAX %ld AVG %ld\n",
(int_t) PS.fill_pelt[4], (int_t)(PS.fill_pelt[3]/PS.fill_pelt[5]));
printf("\t Per elt RL MAX %ld AVG %ld\n",
(int_t) PS.fill_pelt[1], (int_t)(PS.fill_pelt[0]/PS.fill_pelt[2]));
printf("\tM Nops:\t MAX %.2f\tAVG %.2f\n",
mem_glob[11]*1e-6, (stat_glob[16]/nprocs_symb)*1e-6);
printf("\tEXPANSIONS: MAX/AVG\n");
printf("\tTOTAL: %d / %.2f\n",
(int_t) mem_glob[6], stat_glob[11]/nprocs_symb);
printf("\tREALLOC: %.f / %.2f RL_CP %.f / %.2f PR_CP %.f / %.2f\n",
mem_glob[7], stat_glob[12]/nprocs_symb,
mem_glob[8], stat_glob[13]/nprocs_symb,
mem_glob[9], stat_glob[14]/nprocs_symb);
printf ("\n\tDATA MSGS noMsgs*10^3 %.3f/%.3f size (MB) %.3f/%.3f \n",
stat_msgs_g[2]*1e-3, stat_msgs_g[4]/nprocs_symb*1e-3,
stat_msgs_g[0]*1e-6, stat_msgs_g[8] / stat_msgs_g[6]*1e-6);
printf ("\tTOTAL MSGS noMsgs*10^3 %.3f/%.3f size (MB) %.3f/%.3f \n",
stat_msgs_g[3]*1e-3, stat_msgs_g[5]/nprocs_symb*1e-3,
stat_msgs_g[1]*1e-6, stat_msgs_g[9]/stat_msgs_g[7]*1e-6);
#endif
#if ( PROFlevel>=1 )
printf("Distribute matrix time = %8.3f\n", t_symbFact[0]);
printf("Count vertices time = %8.3f\n", t_symbFact[2]);
printf("Symbfact DIST time = %8.3f\n", t_symbFact[1]);
printf("\nLvl\t Time\t Init\t Inter\t Intra\n");
time_lvlsg[0] = 0.;
for (i = 0; i < nlvls; i++) {
for (j = 1; j < 9; j++)
time_lvlsg[j] = 0.;
for (p = 0; p < nprocs_symb; p++) {
k = p * 3 * nlvls;
t = time_lvlsT[i*3+k] + time_lvlsT[i*3+k+1] + time_lvlsT[i*3+k+2];
if (t > time_lvlsg[1]) {
time_lvlsg[1] = t; j = p;
}
time_lvlsg[2] += t;
if (time_lvlsT[i*3+k] > time_lvlsg[3])
time_lvlsg[3] = time_lvlsT[i*3+k];
time_lvlsg[4] += time_lvlsT[i*3+k];
if (time_lvlsT[i*3+k+1] > time_lvlsg[5])
time_lvlsg[5] = time_lvlsT[i*3+k+1];
time_lvlsg[6] += time_lvlsT[i*3+k+1];
if (time_lvlsT[i*3+k+2] > time_lvlsg[7])
time_lvlsg[7] = time_lvlsT[i*3+k+2];
time_lvlsg[8] += time_lvlsT[i*3+k+2];
}
time_lvlsg[0] += time_lvlsg[1];
printf ("%d \t%.3f/%.3f\t%.3f/%.3f\t%.3f/%.3f\t%.3f/%.3f\n", i,
time_lvlsg[1], time_lvlsg[2] / nprocs_symb,
time_lvlsg[3], time_lvlsg[4] / nprocs_symb,
time_lvlsg[5], time_lvlsg[6] /nprocs_symb,
time_lvlsg[7], time_lvlsg[8] / nprocs_symb);
}
printf("\t %8.3f \n", time_lvlsg[0]);
#endif
}
#endif
#if ( PROFlevel>=1 )
SUPERLU_FREE (time_lvls);
SUPERLU_FREE (time_lvlsT);
#endif
symbfact_free (iam, nprocs_symb, &Llu_symbfact, &VInfo, &CS);
} /* if (iam < nprocs_symb) */
else {
/* update Pslu_freeable before returning */
Pslu_freeable->nvtcs_loc = 0;
Pslu_freeable->xlsub = NULL; Pslu_freeable->lsub = NULL;
Pslu_freeable->xusub = NULL; Pslu_freeable->usub = NULL;
Pslu_freeable->supno_loc = NULL;
Pslu_freeable->xsup_beg_loc = NULL;
Pslu_freeable->xsup_end_loc = NULL;
SUPERLU_FREE( tempArray );
PS.allocMem -= n * sizeof(int_t);
}
if (iam < nprocs_symb && nprocs_symb != 1)
freeComm (iam, nprocs_symb, commLvls, symb_comm);
if (commLvls != NULL)
SUPERLU_FREE( commLvls );
#if ( DEBUGlevel>=1 )
CHECK_MALLOC(iam, "Exit psymbfact()");
#endif
return (- PS.allocMem);
} /* SYMBFACT_DIST */
static int_t
initParmsAndStats
(
psymbfact_stat_t *PS /* Output -statistics*/
)
/*! \brief
* <pre>
* Purpose
* =======
* Initialize relaxation parameters and statistics variables
* </pre>
*/
{
int i;
PS->nDnsCurSep = 0;
PS->nDnsUpSeps = 0;
PS->relax_gen = 1.0;
PS->relax_curSep = 1.0;
PS->relax_seps = 1.0;
PS->fill_par = sp_ienv_dist(6);
PS->nops = 0.;
PS->no_shmSnd = 0.;
PS->no_msgsSnd = 0.;
PS->maxsz_msgSnd = 0;
PS->sz_msgsSnd = 0.;
PS->no_shmRcvd = 0.;
PS->no_msgsRcvd = 0.;
PS->maxsz_msgRcvd = 0;
PS->sz_msgsRcvd = 0.;
PS->no_msgsCol = 0.;
PS->maxsz_msgCol = 0;
PS->sz_msgsCol = 0.;
for (i = 0; i < 6; i++)
PS->fill_pelt[i] = 0.;
PS->estimUSz = 0;
PS->estimLSz = 0;
PS->maxSzLPr = 0;
PS->maxSzUPr = 0;
PS->maxSzBuf = 0;
PS->szDnsSep = 0;
PS->allocMem = 0;
return 0;
}
static float
cntsVtcs
(
int_t n, /* Input - order of the input matrix */
int iam, /* Input - my processor number */
int nprocs_symb, /* Input - no of processors for symbolic factorization */
Pslu_freeable_t *Pslu_freeable, /* Input -globToLoc and maxNvtcsPProc */
Llu_symbfact_t *Llu_symbfact, /* Input/Output -local L, U data structures */
vtcsInfo_symbfact_t *VInfo, /* Input - local info on vertices distribution */
int_t *tempArray, /* Input - temporary storage */
int_t *fstVtxSep, /* Input - first vertex of each node in the tree */
int_t *sizes, /* Input - sizes of each node in the tree */
psymbfact_stat_t *PS, /* Input/Output -statistics */
MPI_Comm *commLvls
)
/*! \brief
*
* <pre>
* Purpose
* =======
*
* Computes an estimation of the number of elements in columns of L
* and rows of U. Stores this information in cntelt_vtcs, and it will
* be used in the right-looking symbolic factorization.
* </pre>
*/
{
int fstP, lstP, szSep, npNode, i, j;
int_t nvtcs_loc, ind_blk, vtx, vtx_lid, ii, jj, lv, vtx_elt, cur_blk;
int_t fstVtx, lstVtx, fstVtx_blk, lstVtx_blk;
int_t nelts, nelts_new_blk;
int_t *xlsub, *lsub, *xusub, *usub, *globToLoc, maxNvtcsPProc;
int_t *minElt_vtx, *cntelt_vtcs;
/* Initialization */
xlsub = Llu_symbfact->xlsub; lsub = Llu_symbfact->lsub;
xusub = Llu_symbfact->xusub; usub = Llu_symbfact->usub;
cntelt_vtcs = Llu_symbfact->cntelt_vtcs;
globToLoc = Pslu_freeable->globToLoc;
nvtcs_loc = VInfo->nvtcs_loc;
maxNvtcsPProc = Pslu_freeable->maxNvtcsPProc;
if (Llu_symbfact->szLsub - VInfo->nnz_ainf_loc > n)
minElt_vtx = lsub;
else {
/* allocate memory for minElt_vtx */
if (!(minElt_vtx = intMalloc_dist(n))) {
fprintf(stderr, "Malloc fails for minElt_vtx[].");
return (PS->allocMem);
}
PS->allocMem += n * sizeof (int_t);
}
for (ii = 0; ii < n; ii++)
tempArray[ii] = n;
for (ii = 0; ii < nvtcs_loc; ii++)
cntelt_vtcs[ii] = 0;
szSep = nprocs_symb;
i = 0;
cur_blk = 0;
vtx_lid = 0;
while (szSep >= 1) {
/* for each level in the separator tree */
npNode = nprocs_symb / szSep;
fstP = 0;
/* for each node in the level */
for (j = i; j < i + szSep; j++) {
fstVtx = fstVtxSep[j];
lstVtx = fstVtx + sizes[j];
lstP = fstP + npNode;
if (fstP <= iam && iam < lstP) {
ind_blk = cur_blk;
ii = vtx_lid;
while (VInfo->begEndBlks_loc[ind_blk] < lstVtx &&
ind_blk < 2 * VInfo->nblks_loc) {
fstVtx_blk = VInfo->begEndBlks_loc[ind_blk];
lstVtx_blk = VInfo->begEndBlks_loc[ind_blk + 1];
ind_blk += 2;
for (vtx = fstVtx_blk; vtx < lstVtx_blk; vtx++, ii++) {
for (jj = xlsub[ii]; jj < xlsub[ii+1]; jj++) {
vtx_elt = lsub[jj];
if (tempArray[vtx_elt] == n) {
tempArray[vtx_elt] = vtx;
}
}
for (jj = xusub[ii]; jj < xusub[ii+1]; jj++) {
vtx_elt = usub[jj];
if (tempArray[vtx_elt] == n) {
tempArray[vtx_elt] = vtx;
}
}
}
}
if (szSep == nprocs_symb)
vtx_lid = ii;
else {
MPI_Allreduce (&(tempArray[fstVtx]), &(minElt_vtx[fstVtx]),
(int) (n - fstVtx), mpi_int_t, MPI_MIN, commLvls[j]);
#if ( PRNTlevel>=1 )
PS->no_msgsCol += (float) (2 * (int) LOG2( npNode ));
PS->sz_msgsCol += (float) (n - fstVtx);
if (PS->maxsz_msgCol < n - fstVtx)
PS->maxsz_msgCol = n - fstVtx;
#endif
nelts = 0;
for (ii = fstVtx; ii < lstVtx; ii++)
tempArray[ii] = 0;
for (ii = fstVtx; ii < n; ii++) {
if (minElt_vtx[ii] != n) {
if (minElt_vtx[ii] < fstVtx)
nelts ++;
else
tempArray[minElt_vtx[ii]] ++;
if (ii > lstVtx)
tempArray[ii] = minElt_vtx[ii];
}
}
ind_blk = cur_blk;
lv = fstVtx;
while (VInfo->begEndBlks_loc[ind_blk] < lstVtx &&
ind_blk < 2 * VInfo->nblks_loc) {
fstVtx_blk = VInfo->begEndBlks_loc[ind_blk];
lstVtx_blk = VInfo->begEndBlks_loc[ind_blk + 1];
ind_blk += 2;
for (ii = lv; ii < fstVtx_blk; ii++)
nelts += tempArray[ii];
lv = lstVtx_blk;
nelts_new_blk = 0;
for (vtx = fstVtx_blk; vtx < lstVtx_blk; vtx++, vtx_lid++) {
nelts_new_blk += tempArray[vtx];
cntelt_vtcs[vtx_lid] = nelts;
}
nelts += nelts_new_blk;
}
} /* if (szSep != nprocs_symb) */
cur_blk = ind_blk;
}
fstP += npNode;
}
i += szSep;
szSep = szSep / 2;
}
/* free memory */
if (minElt_vtx != lsub) {
SUPERLU_FREE (minElt_vtx);
PS->allocMem -= n * sizeof(int_t);
}
return (SUCCES_RET);
}
static float
symbfact_mapVtcs
(
int iam, /* Input -process number */
int nprocs_num, /* Input -number of processors */
int nprocs_symb, /* Input -number of procs for symbolic factorization */
SuperMatrix *A, /* Input -input distributed matrix A */
int_t *fstVtxSep, /* Input -first vertex in each separator */
int_t *sizes, /* Input -size of each separator in the separator tree */
Pslu_freeable_t *Pslu_freeable, /* Output -globToLoc and maxNvtcsPProc
computed */
vtcsInfo_symbfact_t *VInfo, /* Output -local info on vertices distribution */
int_t *tempArray, /* Input -temp array of size n = order of the matrix */
int_t maxSzBlk, /* Input -maximum number of vertices in a block */
psymbfact_stat_t *PS /* Input/Output -statistics */
)
{
/*! \brief
*
* <pre>
* Purpose
* =======
*
* symbfact_mapVtcs maps the vertices of the graph of the input
* matrix A on nprocs_symb processors, using the separator tree
* returned by a graph partitioning algorithm from the previous step
* of the symbolic factorization. The number of processors
* nprocs_symb must be a power of 2.
*