Actual source code: superlu_dist.c
1: /*$Id: superlu_DIST.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2: /*
3: Provides an interface to the SuperLU_DIST sparse solver
4: Usage:
5: mpirun -np <procs> main -mat_aij_superlu_dist -r <proc rows> -c <proc columns>
6: or
7: mpirun -np <procs> main -mat_aij_superlu_dist (use the default process grid)
9: Command line options:
10: -mat_aij_superlu_dist_equil <YES/NO>
11: -mat_aij_superlu_dist_rowperm <NATURAL/LargeDiag>
12: -mat_aij_superlu_dist_colperm <NATURAL/COLAMD/MMD_ATA/MMD_AT_PLUS_A>
13: -mat_aij_superlu_dist_replacetinypivot <YES/NO>
14: -mat_aij_superlu_dist_iterrefine <NO/DOUBLE>
15: -mat_aij_superlu_dist_statprint <YES/NO>
17: SuperLU_DIST default options:
18: equil: YES
19: rowperm: LargeDiag
20: colperm: MMD_AT_PLUS_A
21: replacetinypivot: YES
22: iterrefine: NO
23: statprint: YES
24: */
26: #include src/mat/impls/aij/seq/aij.h
27: #include src/mat/impls/aij/mpi/mpiaij.h
29: #if defined(PETSC_HAVE_SUPERLUDIST) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
31: EXTERN_C_BEGIN
32: #include "superlu_ddefs.h"
33: EXTERN_C_END
36: typedef struct {
37: gridinfo_t grid;
38: superlu_options_t options;
39: SuperMatrix A_sup;
40: ScalePermstruct_t ScalePermstruct;
41: LUstruct_t LUstruct;
42: int StatPrint;
43: } Mat_MPIAIJ_SuperLU_DIST;
45: extern int MatDestroy_MPIAIJ(Mat);
46: extern int MatDestroy_SeqAIJ(Mat);
48: #if !defined(PETSC_HAVE_SUPERLU)
49: /* SuperLU function: Convert a row compressed storage into a column compressed storage. */
50: void dCompRow_to_CompCol(int m, int n, int nnz,
51: double *a, int *colind, int *rowptr,
52: double **at, int **rowind, int **colptr)
53: {
54: register int i, j, col, relpos;
55: int *marker;
57: /* Allocate storage for another copy of the matrix. */
58: *at = (double *) doubleMalloc_dist(nnz);
59: *rowind = (int *) intMalloc_dist(nnz);
60: *colptr = (int *) intMalloc_dist(n+1);
61: marker = (int *) intCalloc_dist(n);
63: /* Get counts of each column of A, and set up column pointers */
64: for (i = 0; i < m; ++i)
65: for (j = rowptr[i]; j < rowptr[i+1]; ++j) ++marker[colind[j]];
66: (*colptr)[0] = 0;
67: for (j = 0; j < n; ++j) {
68: (*colptr)[j+1] = (*colptr)[j] + marker[j];
69: marker[j] = (*colptr)[j];
70: }
72: /* Transfer the matrix into the compressed column storage. */
73: for (i = 0; i < m; ++i) {
74: for (j = rowptr[i]; j < rowptr[i+1]; ++j) {
75: col = colind[j];
76: relpos = marker[col];
77: (*rowind)[relpos] = i;
78: (*at)[relpos] = a[j];
79: ++marker[col];
80: }
81: }
83: SUPERLU_FREE(marker);
84: }
85: #else
86: EXTERN_C_BEGIN
87: extern void dCompRow_to_CompCol(int,int,int,double*,int*,int*,double**,int**,int**);
88: EXTERN_C_END
89: #endif /* PETSC_HAVE_SUPERLU*/
91: extern int MatDestroy_MPIAIJ_SuperLU_DIST(Mat A)
92: {
93: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
94: Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)a->spptr;
95: int ierr, size=a->size;
96:
98: /* Deallocate SuperLU_DIST storage */
99: Destroy_CompCol_Matrix_dist(&lu->A_sup);
100: Destroy_LU(A->N, &lu->grid, &lu->LUstruct);
101: ScalePermstructFree(&lu->ScalePermstruct);
102: LUstructFree(&lu->LUstruct);
104: /* Release the SuperLU_DIST process grid. */
105: superlu_gridexit(&lu->grid);
107: PetscFree(lu);
108:
109: if (size == 1){
110: MatDestroy_SeqAIJ(A);
111: } else {
112: MatDestroy_MPIAIJ(A);
113: }
114:
115: return(0);
116: }
118: extern int MatSolve_MPIAIJ_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
119: {
120: Mat_MPIAIJ *aa = (Mat_MPIAIJ*)A->data;
121: Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)aa->spptr;
122: int ierr, size=aa->size;
123: int_t m=A->M, N=A->N;
124: superlu_options_t options=lu->options;
125: SuperLUStat_t stat;
126: double berr[1],*bptr;
127: int info, nrhs=1;
128: Vec x_seq;
129: IS iden;
130: VecScatter scat;
133: if (size > 1) { /* convert mpi vector b to seq vector x_seq */
134: VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);
135: ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);
136: VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);
137: ISDestroy(iden);
139: VecScatterBegin(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);
140: VecScatterEnd(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);
141: VecGetArray(x_seq,&bptr);
142: } else {
143: VecCopy(b_mpi,x);
144: VecGetArray(x,&bptr);
145: }
146:
147: options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only.*/
148:
149: PStatInit(&stat); /* Initialize the statistics variables. */
150: pdgssvx_ABglobal(&options, &lu->A_sup, &lu->ScalePermstruct, bptr, m, nrhs,
151: &lu->grid, &lu->LUstruct, berr, &stat, &info);
152: if (lu->StatPrint) PStatPrint(&stat, &lu->grid); /* Print the statistics. */
153: PStatFree(&stat);
154:
155: if (size > 1) { /* convert seq x to mpi x */
156: VecRestoreArray(x_seq,&bptr);
157: VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);
158: VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);
159: VecScatterDestroy(scat);
160: VecDestroy(x_seq);
161: } else {
162: VecRestoreArray(x,&bptr);
163: }
165: return(0);
166: }
168: extern int MatLUFactorNumeric_MPIAIJ_SuperLU_DIST(Mat A,Mat *F)
169: {
170: Mat_MPIAIJ *fac = (Mat_MPIAIJ*)(*F)->data;
171: Mat *tseq,A_seq = PETSC_NULL;
172: Mat_SeqAIJ *aa;
173: Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)fac->spptr;
174: int M=A->M,N=A->N,info,ierr,size=fac->size;
175: SuperLUStat_t stat;
176: double *berr=0, *bptr=0;
177: int_t *asub, *xa;
178: double *a;
179: SuperMatrix A_sup;
180: IS isrow,iscol;
183: if (size > 1) { /* convert mpi A to seq mat A */
184: ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);
185: ISCreateStride(PETSC_COMM_SELF,N,0,1,&iscol);
187: MatGetSubMatrices(A,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&tseq);
189: ISDestroy(isrow);
190: ISDestroy(iscol);
191: A_seq = *tseq;
192: PetscFree(tseq);
193:
194: aa = (Mat_SeqAIJ*)A_seq->data;
195: } else {
196: aa = (Mat_SeqAIJ*)A->data;
197: }
199: /* Allocate storage for compressed column representation. */
200: dallocateA_dist(N, aa->nz, &a, &asub, &xa);
201:
202: /* Convert Petsc NR matrix storage to SuperLU_DIST NC storage */
203: dCompRow_to_CompCol(M,N,aa->nz,aa->a,aa->j,aa->i,&a, &asub, &xa);
205: /* Create compressed column matrix A_sup. */
206: dCreate_CompCol_Matrix_dist(&A_sup, M, N, aa->nz, a, asub, xa, NC, D, GE);
208: /* Factor the matrix. */
209: PStatInit(&stat); /* Initialize the statistics variables. */
210: pdgssvx_ABglobal(&lu->options, &A_sup, &lu->ScalePermstruct, bptr, M, 0,
211: &lu->grid, &lu->LUstruct, berr, &stat, &info);
212: if (lu->StatPrint) PStatPrint(&stat, &lu->grid); /* Print the statistics. */
213: PStatFree(&stat);
215: lu->A_sup = A_sup;
216: lu->options.Fact = SamePattern; /* Sparsity pattern of A and perm_c can be reused. */
217: if (size > 1){
218: MatDestroy(A_seq);
219: }
220:
221: return(0);
222: }
224: /* Note the Petsc r and c permutations are ignored */
226: extern int MatLUFactorSymbolic_MPIAIJ_SuperLU_DIST(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
227: {
228: Mat_MPIAIJ *fac;
229: Mat_MPIAIJ_SuperLU_DIST *lu; /* ptr to Mat_MPIAIJ_SuperLU_DIST */
230: int ierr,M=A->M,N=A->N,size;
231: int_t nprow, npcol;
232: gridinfo_t grid;
233: superlu_options_t options;
234: ScalePermstruct_t ScalePermstruct;
235: LUstruct_t LUstruct;
236: char opt[256];
237: PetscTruth flg,flg1;
240: /* Initialize the SuperLU process grid. */
241: MPI_Comm_size(PETSC_COMM_WORLD,&size);
242: nprow = size/2; /* Default process rows. */
243: if (nprow == 0) nprow = 1;
244: npcol = size/nprow; /* Default process columns. */
245:
246: PetscOptionsGetInt(PETSC_NULL,"-r",&nprow,PETSC_NULL);
247: PetscOptionsGetInt(PETSC_NULL,"-c",&npcol,PETSC_NULL);
248:
249: if ( size != nprow * npcol )
250: SETERRQ(1,"Number of processes should be equal to nprow*npcol");
251:
252: superlu_gridinit(MPI_COMM_WORLD, nprow, npcol, &grid);
254: /* Create the factorization matrix F */
255: MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,M,N,0,PETSC_NULL,0,PETSC_NULL,F);
256:
258: (*F)->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_SuperLU_DIST;
259: (*F)->ops->solve = MatSolve_MPIAIJ_SuperLU_DIST;
260: (*F)->ops->destroy = MatDestroy_MPIAIJ_SuperLU_DIST;
261: (*F)->factor = FACTOR_LU;
262: fac = (Mat_MPIAIJ*)(*F)->data;
264: ierr = PetscNew(Mat_MPIAIJ_SuperLU_DIST,&lu);
265: fac->spptr = (void*)lu;
267: /* Set the input options */
268: set_default_options(&options);
269: options.IterRefine = NOREFINE;
271: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_equil",opt,256,&flg);
272: if (flg) {
273: PetscStrcmp(opt,"NO",&flg1);
274: if (flg1) options.Equil = NO;
275: }
277: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_rowperm",opt,256,&flg);
278: if (flg) {
279: PetscStrcmp(opt,"NATURAL",&flg1);
280: if (flg1) options.RowPerm = NOROWPERM;
281: }
283: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_colperm",opt,256,&flg);
284: while (flg) {
285: PetscStrcmp(opt,"NATURAL",&flg1);
286: if (flg1) {
287: options.ColPerm = NATURAL;
288: break;
289: }
290: PetscStrcmp(opt,"MMD_ATA",&flg1);
291: if (flg1) {
292: options.ColPerm = MMD_ATA;
293: break;
294: }
295: PetscStrcmp(opt,"COLAMD",&flg1);
296: if (flg1) {
297: options.ColPerm = COLAMD;
298: break;
299: }
300: break;
301: }
303: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_replacetinypivot",opt,256,&flg);
304: if (flg) {
305: PetscStrcmp(opt,"NO",&flg1);
306: if (flg1) options.ReplaceTinyPivot = NO;
307: }
309: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_iterrefine",opt,256,&flg);
310: if (flg) {
311: PetscStrcmp(opt,"DOUBLE",&flg1);
312: if (flg1) options.IterRefine = DOUBLE;
313: }
315: lu->StatPrint = 1;
316: PetscOptionsGetString(PETSC_NULL,"-mat_aij_superlu_dist_statprint",opt,256,&flg);
317: if (flg) {
318: PetscStrcmp(opt,"NO",&flg1);
319: if (flg1) {
320: lu->StatPrint = 0;
321: }
322: }
324: /* Initialize ScalePermstruct and LUstruct. */
325: ScalePermstructInit(M, N, &ScalePermstruct);
326: LUstructInit(M, N, &LUstruct);
328: lu->ScalePermstruct = ScalePermstruct;
329: lu->LUstruct = LUstruct;
330: lu->options = options;
331: lu->grid = grid;
332: fac->size = size;
334: return(0);
335: }
337: int MatUseSuperLU_DIST_MPIAIJ(Mat A)
338: {
340: A->ops->lufactorsymbolic = MatLUFactorSymbolic_MPIAIJ_SuperLU_DIST;
341: A->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_SuperLU_DIST;
342: return(0);
343: }
345: #else
347: int MatUseSuperLU_DIST_MPIAIJ(Mat A)
348: {
350: return(0);
351: }
353: #endif