forked from xiaoyeli/superlu_dist
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpsdrive3d1.c
executable file
·464 lines (412 loc) · 13.8 KB
/
psdrive3d1.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
/*! \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 Driver program for PSGSSVX3D example
*
* <pre>
* -- Distributed SuperLU routine (version 9.0) --
* Lawrence Berkeley National Lab, Georgia Institute of Technology,
* Oak Ridge National Lab
* September 10, 2021
*
*/
#include "superlu_sdefs.h"
/*! \brief
*
* <pre>
* Purpose
* =======
*
* The driver program PSDRIVE3D1.
*
* This example illustrates how to use PSGSSVX3D to sovle the systems
* with the same A but different right-hand side, possibly with
* different number of right-hand sides.
* In this case, we factorize A only once in the first call to PSGSSVX3D,
* and reuse the following data structures in the subsequent call to
* PSGSSVX3D:
* ScalePermstruct : DiagScale, R, C, perm_r, perm_c
* LUstruct : Glu_persist, Llu
* SOLVEstruct : communication metadata for SpTRSV, SpMV, and
* 3D<->2D gather/scatter of {A,B} stored in A3d.
*
* The program may be run by typing:
* mpiexec -np <p> psdrive3d1 -r <proc rows> -c <proc columns> \
* -d <proc Z-dimension> <input_file>
* NOTE: total number of processes p = r * c * d
* (d must be a power-of-two, e.g., 1, 2, 4, ...)
*
* </pre>
*/
static void matCheck(int n, int m, float* A, int LDA,
float* B, int LDB)
{
for(int j=0; j<m;j++)
for (int i = 0; i < n; ++i) {
assert(A[i+ LDA*j] == B[i+ LDB*j]);
}
printf("B check passed\n");
return;
}
static void checkNRFMT(NRformat_loc*A, NRformat_loc*B)
{
/*
int_t nnz_loc;
int_t m_loc;
int_t fst_row;
void *nzval;
int_t *rowptr;
int_t *colind;
*/
assert(A->nnz_loc == B->nnz_loc);
assert(A->m_loc == B->m_loc);
assert(A->fst_row == B->fst_row);
#if 0
double *Aval = (double *)A->nzval, *Bval = (double *)B->nzval;
Printdouble5("A", A->nnz_loc, Aval);
Printdouble5("B", B->nnz_loc, Bval);
fflush(stdout);
#endif
float * Aval = (float *) A->nzval;
float * Bval = (float *) B->nzval;
for (int_t i = 0; i < A->nnz_loc; i++)
{
assert( Aval[i] == Bval[i] );
assert((A->colind)[i] == (B->colind)[i]);
printf("colind[] correct\n");
}
for (int_t i = 0; i < A->m_loc + 1; i++)
{
assert((A->rowptr)[i] == (B->rowptr)[i]);
}
printf("Matrix check passed\n");
}
int
main (int argc, char *argv[])
{
superlu_dist_options_t options;
SuperLUStat_t stat;
SuperMatrix A; // Now, A is on all 3D processes
sScalePermstruct_t ScalePermstruct;
sLUstruct_t LUstruct;
sSOLVEstruct_t SOLVEstruct;
gridinfo3d_t grid;
float *berr;
float *b, *xtrue, *b1, *b2;
int m, n, i, j, m_loc;
int nprow, npcol, npdep;
int lookahead, colperm, rowperm, ir;
int iam, info, ldb, ldx, nrhs;
char **cpp, c, *suffix;
FILE *fp, *fopen ();
extern int cpp_defs ();
int ii, omp_mpi_level;
nprow = 1; /* Default process rows. */
npcol = 1; /* Default process columns. */
npdep = 1; /* replication factor must be power of two */
nrhs = 1; /* Number of right-hand side. */
lookahead = -1;
colperm = -1;
rowperm = -1;
ir = -1;
/* ------------------------------------------------------------
INITIALIZE MPI ENVIRONMENT.
------------------------------------------------------------ */
// MPI_Init (&argc, &argv);
int required = MPI_THREAD_MULTIPLE;
int provided;
MPI_Init_thread(&argc, &argv, required, &provided);
if (provided < required)
{
int rank;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
if (!rank) {
printf("The MPI library doesn't provide MPI_THREAD_MULTIPLE \n");
printf("\tprovided omp_mpi_level: %d\n", provided);
}
}
/* Parse command line argv[]. */
for (cpp = argv + 1; *cpp; ++cpp)
{
if (**cpp == '-')
{
c = *(*cpp + 1);
++cpp;
switch (c)
{
case 'h':
printf ("Options:\n");
printf ("\t-r <int>: process rows (default %d)\n", nprow);
printf ("\t-c <int>: process columns (default %d)\n", npcol);
printf ("\t-d <int>: process Z-dimension (default %d)\n", npdep);
exit (0);
break;
case 'r':
nprow = atoi (*cpp);
break;
case 'c':
npcol = atoi (*cpp);
break;
case 'd':
npdep = atoi (*cpp);
break;
case 'l': lookahead = atoi(*cpp);
break;
case 'p': rowperm = atoi(*cpp);
break;
case 'q': colperm = atoi(*cpp);
break;
case 'i': ir = atoi(*cpp);
break;
}
}
else
{ /* Last arg is considered a filename */
if (!(fp = fopen (*cpp, "r")))
{
ABORT ("File does not exist");
}
break;
}
}
/* ------------------------------------------------------------
INITIALIZE THE SUPERLU PROCESS GRID.
------------------------------------------------------------ */
superlu_gridinit3d (MPI_COMM_WORLD, nprow, npcol, npdep, &grid);
if(grid.iam==0) {
MPI_Query_thread(&omp_mpi_level);
switch (omp_mpi_level) {
case MPI_THREAD_SINGLE:
printf("MPI_Query_thread with MPI_THREAD_SINGLE\n");
fflush(stdout);
break;
case MPI_THREAD_FUNNELED:
printf("MPI_Query_thread with MPI_THREAD_FUNNELED\n");
fflush(stdout);
break;
case MPI_THREAD_SERIALIZED:
printf("MPI_Query_thread with MPI_THREAD_SERIALIZED\n");
fflush(stdout);
break;
case MPI_THREAD_MULTIPLE:
printf("MPI_Query_thread with MPI_THREAD_MULTIPLE\n");
fflush(stdout);
break;
}
fflush(stdout);
}
/* Bail out if I do not belong in the grid. */
iam = grid.iam;
if (iam == -1) goto out;
if (!iam) {
int v_major, v_minor, v_bugfix;
#ifdef __INTEL_COMPILER
printf("__INTEL_COMPILER is defined\n");
#endif
printf("__STDC_VERSION__ %ld\n", __STDC_VERSION__);
superlu_dist_GetVersionNumber(&v_major, &v_minor, &v_bugfix);
printf("Library version:\t%d.%d.%d\n", v_major, v_minor, v_bugfix);
printf("Input matrix file:\t%s\n", *cpp);
printf("3D process grid: %d X %d X %d\n", nprow, npcol, npdep);
//printf("2D Process grid: %d X %d\n", (int)grid.nprow, (int)grid.npcol);
fflush(stdout);
}
#if ( DEBUGlevel>=1 )
CHECK_MALLOC (iam, "Enter main()");
#endif
/* ------------------------------------------------------------
GET THE MATRIX FROM FILE AND SETUP THE RIGHT HAND SIDE.
------------------------------------------------------------ */
for (ii = 0; ii<strlen(*cpp); ii++) {
if((*cpp)[ii]=='.'){
suffix = &((*cpp)[ii+1]);
// printf("%s\n", suffix);
}
}
// *fp0 = *fp;
screate_matrix_postfix3d(&A, nrhs, &b, &ldb,
&xtrue, &ldx, fp, suffix, &(grid));
//printf("ldx %d, ldb %d\n", ldx, ldb);
#if 0 // following code is only for checking *Gather* routine
NRformat_loc *Astore, *Astore0;
float* B2d;
NRformat_loc Atmp = dGatherNRformat_loc(
(NRformat_loc *) A.Store,
b, ldb, nrhs, &B2d,
&grid);
Astore = &Atmp;
SuperMatrix Aref;
float *bref, *xtrueref;
if ( grid.zscp.Iam == 0 ) // only in process layer 0
{
screate_matrix_postfix(&Aref, nrhs, &bref, &ldb,
&xtrueref, &ldx, fp0,
suffix, &(grid.grid2d));
Astore0 = (NRformat_loc *) Aref.Store;
/*
if ( (grid.grid2d).iam == 0 ) {
printf(" iam %d\n", 0);
checkNRFMT(Astore, Astore0);
} else if ((grid.grid2d).iam == 1 ) {
printf(" iam %d\n", 1);
checkNRFMT(Astore, Astore0);
}
*/
// bref, xtrueref are created on 2D
matCheck(Astore->m_loc, nrhs, B2d, Astore->m_loc, bref, ldb);
}
// MPI_Finalize(); exit(0);
#endif
/* Save two copies of the RHS */
if ( !(b1 = floatMalloc_dist(ldb * nrhs)) )
ABORT("Malloc fails for b1[]");
if ( !(b2 = floatMalloc_dist(ldb * nrhs)) )
ABORT("Malloc fails for b1[]");
for (j = 0; j < nrhs; ++j) {
for (i = 0; i < ldb; ++i) {
b1[i+j*ldb] = b[i+j*ldb];
b2[i+j*ldb] = b[i+j*ldb];
}
}
if (!(berr = floatMalloc_dist (nrhs)))
ABORT ("Malloc fails for berr[].");
/* ------------------------------------------------------------
1. SOLVE THE LINEAR SYSTEM FOR THE FIRST TIME, WITH 1 RHS.
------------------------------------------------------------*/
/* Set the default input options:
options.Fact = DOFACT;
options.Equil = YES;
options.ParSymbFact = NO;
options.ColPerm = METIS_AT_PLUS_A;
options.RowPerm = LargeDiag_MC64;
options.ReplaceTinyPivot = NO;
options.IterRefine = SLU_DOUBLE;
options.Trans = NOTRANS;
options.SolveInitialized = NO;
options.RefineInitialized = NO;
options.PrintStat = YES;
options->num_lookaheads = 10;
options->lookahead_etree = NO;
options->SymPattern = NO;
options.DiagInv = NO;
*/
set_default_options_dist (&options);
options.Algo3d = YES;
options.IterRefine = SLU_SINGLE;
#if 0
options.RowPerm = NOROWPERM;
options.IterRefine = NOREFINE;
options.ColPerm = NATURAL;
options.Equil = NO;
options.ReplaceTinyPivot = YES;
#endif
if (rowperm != -1) options.RowPerm = rowperm;
if (colperm != -1) options.ColPerm = colperm;
if (lookahead != -1) options.num_lookaheads = lookahead;
if (ir != -1) options.IterRefine = ir;
if (!iam) {
print_options_dist(&options);
fflush(stdout);
}
// matrix is on 3D process grid
m = A.nrow;
n = A.ncol;
/* Initialize ScalePermstruct and LUstruct. */
sScalePermstructInit (m, n, &ScalePermstruct);
sLUstructInit (n, &LUstruct);
/* Initialize the statistics variables. */
PStatInit (&stat);
/* Call the linear equation solver. */
psgssvx3d (&options, &A, &ScalePermstruct, b, ldb, nrhs, &grid,
&LUstruct, &SOLVEstruct, berr, &stat, &info);
if ( info ) { /* Something is wrong */
if ( iam==0 ) {
printf("ERROR: INFO = %d returned from psgssvx3d()\n", info);
fflush(stdout);
}
} else {
/* Check the accuracy of the solution. */
if ( !iam ) printf("\tSolve the first system:\n");
psinf_norm_error (iam, ((NRformat_loc *) A.Store)->m_loc,
nrhs, b, ldb, xtrue, ldx, grid.comm);
}
if ( grid.zscp.Iam == 0 ) { // process layer 0
PStatPrint (&options, &stat, &(grid.grid2d)); /* Print 2D statistics.*/
}
PStatFree (&stat);
fflush(stdout);
/* ------------------------------------------------------------
2. NOW SOLVE ANOTHER SYSTEM WITH THE SAME A BUT DIFFERENT
RIGHT-HAND SIDE, WE WILL USE THE EXISTING L AND U FACTORS IN
LUSTRUCT OBTAINED FROM A PREVIOUS FATORIZATION.
------------------------------------------------------------*/
options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
PStatInit(&stat); /* Initialize the statistics variables. */
nrhs = 1;
psgssvx3d (&options, &A, &ScalePermstruct, b1, ldb, nrhs, &grid,
&LUstruct, &SOLVEstruct, berr, &stat, &info);
if ( info ) { /* Something is wrong */
if ( iam==0 ) {
printf("ERROR: INFO = %d returned from psgssvx3d()\n", info);
fflush(stdout);
}
} else {
/* Check the accuracy of the solution. */
if ( !iam ) printf("\tSolve the system with a different B:\n");
psinf_norm_error (iam, ((NRformat_loc *) A.Store)->m_loc,
nrhs, b1, ldb, xtrue, ldx, grid.comm);
}
/* ------------------------------------------------------------
DEALLOCATE STORAGE.
------------------------------------------------------------ */
if ( grid.zscp.Iam == 0 ) { // process layer 0
PStatPrint (&options, &stat, &(grid.grid2d)); /* Print 2D statistics.*/
}
sSolveFinalize (&options, &SOLVEstruct);
sDestroy_LU (n, &(grid.grid2d), &LUstruct);
sDestroy_A3d_gathered_on_2d(&SOLVEstruct, &grid);
Destroy_CompRowLoc_Matrix_dist (&A);
SUPERLU_FREE (b);
SUPERLU_FREE (b1);
SUPERLU_FREE (b2);
SUPERLU_FREE (xtrue);
SUPERLU_FREE (berr);
sScalePermstructFree (&ScalePermstruct);
sLUstructFree (&LUstruct);
PStatFree (&stat);
fclose(fp);
/* ------------------------------------------------------------
RELEASE THE SUPERLU PROCESS GRID.
------------------------------------------------------------ */
out:
superlu_gridexit3d (&grid);
/* ------------------------------------------------------------
TERMINATES THE MPI EXECUTION ENVIRONMENT.
------------------------------------------------------------ */
MPI_Finalize ();
#if ( DEBUGlevel>=1 )
CHECK_MALLOC (iam, "Exit main()");
#endif
}
int
cpp_defs ()
{
printf (".. CPP definitions:\n");
#if ( PRNTlevel>=1 )
printf ("\tPRNTlevel = %d\n", PRNTlevel);
#endif
#if ( DEBUGlevel>=1 )
printf ("\tDEBUGlevel = %d\n", DEBUGlevel);
#endif
#if ( PROFlevel>=1 )
printf ("\tPROFlevel = %d\n", PROFlevel);
#endif
printf ("....\n");
return 0;
}