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_2.0 sparse solver
4: */
6: #include "src/mat/impls/aij/seq/aij.h"
7: #include "src/mat/impls/aij/mpi/mpiaij.h"
8: #if defined(PETSC_HAVE_STDLIB_H) /* This is to get arround weird problem with SuperLU on cray */
9: #include "stdlib.h"
10: #endif
12: EXTERN_C_BEGIN
13: #if defined(PETSC_USE_COMPLEX)
14: #include "superlu_zdefs.h"
15: #else
16: #include "superlu_ddefs.h"
17: #endif
18: EXTERN_C_END
20: typedef enum { GLOBAL,DISTRIBUTED
21: } SuperLU_MatInputMode;
23: typedef struct {
24: int_t nprow,npcol,*row,*col;
25: gridinfo_t grid;
26: superlu_options_t options;
27: SuperMatrix A_sup;
28: ScalePermstruct_t ScalePermstruct;
29: LUstruct_t LUstruct;
30: int StatPrint;
31: int MatInputMode;
32: SOLVEstruct_t SOLVEstruct;
33: MatStructure flg;
34: MPI_Comm comm_superlu;
35: #if defined(PETSC_USE_COMPLEX)
36: doublecomplex *val;
37: #else
38: double *val;
39: #endif
41: /* A few function pointers for inheritance */
42: int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
43: int (*MatView)(Mat,PetscViewer);
44: int (*MatAssemblyEnd)(Mat,MatAssemblyType);
45: int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
46: int (*MatDestroy)(Mat);
48: /* Flag to clean up (non-global) SuperLU objects during Destroy */
49: PetscTruth CleanUpSuperLU_Dist;
50: } Mat_SuperLU_DIST;
52: EXTERN int MatDuplicate_SuperLU_DIST(Mat,MatDuplicateOption,Mat*);
54: EXTERN_C_BEGIN
57: int MatConvert_SuperLU_DIST_Base(Mat A,const MatType type,Mat *newmat) {
58: int ierr;
59: Mat B=*newmat;
60: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST *)A->spptr;
63: if (B != A) {
64: MatDuplicate(A,MAT_COPY_VALUES,&B);
65: }
66: /* Reset the original function pointers */
67: B->ops->duplicate = lu->MatDuplicate;
68: B->ops->view = lu->MatView;
69: B->ops->assemblyend = lu->MatAssemblyEnd;
70: B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
71: B->ops->destroy = lu->MatDestroy;
73: PetscObjectChangeTypeName((PetscObject)B,type);
74: PetscFree(lu);
76: *newmat = B;
77: return(0);
78: }
79: EXTERN_C_END
83: int MatDestroy_SuperLU_DIST(Mat A)
84: {
85: int ierr,size;
86: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->spptr;
87:
89: if (lu->CleanUpSuperLU_Dist) {
90: /* Deallocate SuperLU_DIST storage */
91: if (lu->MatInputMode == GLOBAL) {
92: Destroy_CompCol_Matrix_dist(&lu->A_sup);
93: } else {
94: Destroy_CompRowLoc_Matrix_dist(&lu->A_sup);
95: if ( lu->options.SolveInitialized ) {
96: #if defined(PETSC_USE_COMPLEX)
97: zSolveFinalize(&lu->options, &lu->SOLVEstruct);
98: #else
99: dSolveFinalize(&lu->options, &lu->SOLVEstruct);
100: #endif
101: }
102: }
103: Destroy_LU(A->N, &lu->grid, &lu->LUstruct);
104: ScalePermstructFree(&lu->ScalePermstruct);
105: LUstructFree(&lu->LUstruct);
107: /* Release the SuperLU_DIST process grid. */
108: superlu_gridexit(&lu->grid);
109:
110: MPI_Comm_free(&(lu->comm_superlu));
111: }
113: MPI_Comm_size(A->comm,&size);
114: if (size == 1) {
115: MatConvert_SuperLU_DIST_Base(A,MATSEQAIJ,&A);
116: } else {
117: MatConvert_SuperLU_DIST_Base(A,MATMPIAIJ,&A);
118: }
119: (*A->ops->destroy)(A);
120:
121: return(0);
122: }
126: int MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
127: {
128: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->spptr;
129: int ierr, size;
130: int m=A->M, N=A->N;
131: SuperLUStat_t stat;
132: double berr[1];
133: PetscScalar *bptr;
134: int info, nrhs=1;
135: Vec x_seq;
136: IS iden;
137: VecScatter scat;
138:
140: MPI_Comm_size(A->comm,&size);
141: if (size > 1) {
142: if (lu->MatInputMode == GLOBAL) { /* global mat input, convert b to x_seq */
143: VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);
144: ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);
145: VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);
146: ISDestroy(iden);
148: VecScatterBegin(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);
149: VecScatterEnd(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);
150: VecGetArray(x_seq,&bptr);
151: } else { /* distributed mat input */
152: VecCopy(b_mpi,x);
153: VecGetArray(x,&bptr);
154: }
155: } else { /* size == 1 */
156: VecCopy(b_mpi,x);
157: VecGetArray(x,&bptr);
158: }
159:
160: lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only.*/
162: PStatInit(&stat); /* Initialize the statistics variables. */
164: if (lu->MatInputMode == GLOBAL) {
165: #if defined(PETSC_USE_COMPLEX)
166: pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, m, nrhs,
167: &lu->grid, &lu->LUstruct, berr, &stat, &info);
168: #else
169: pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, m, nrhs,
170: &lu->grid, &lu->LUstruct, berr, &stat, &info);
171: #endif
172: } else { /* distributed mat input */
173: #if defined(PETSC_USE_COMPLEX)
174: pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, (doublecomplex*)bptr, A->M, nrhs, &lu->grid,
175: &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
176: if (info) SETERRQ1(1,"pzgssvx fails, info: %d\n",info);
177: #else
178: pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, bptr, A->M, nrhs, &lu->grid,
179: &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
180: if (info) SETERRQ1(1,"pdgssvx fails, info: %d\n",info);
181: #endif
182: }
183: if (lu->StatPrint) {
184: PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */
185: }
186: PStatFree(&stat);
187:
188: if (size > 1) {
189: if (lu->MatInputMode == GLOBAL){ /* convert seq x to mpi x */
190: VecRestoreArray(x_seq,&bptr);
191: VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);
192: VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);
193: VecScatterDestroy(scat);
194: VecDestroy(x_seq);
195: } else {
196: VecRestoreArray(x,&bptr);
197: }
198: } else {
199: VecRestoreArray(x,&bptr);
200: }
201: return(0);
202: }
206: int MatLUFactorNumeric_SuperLU_DIST(Mat A,Mat *F)
207: {
208: Mat *tseq,A_seq = PETSC_NULL;
209: Mat_SeqAIJ *aa,*bb;
210: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)(*F)->spptr;
211: int M=A->M,N=A->N,info,ierr,size,rank,i,*ai,*aj,*bi,*bj,nz,rstart,*garray,
212: m=A->m, irow,colA_start,j,jcol,jB,countA,countB,*bjj,*ajj;
213: SuperLUStat_t stat;
214: double *berr=0;
215: IS isrow;
216: PetscLogDouble time0,time,time_min,time_max;
217: #if defined(PETSC_USE_COMPLEX)
218: doublecomplex *av, *bv;
219: #else
220: double *av, *bv;
221: #endif
224: MPI_Comm_size(A->comm,&size);
225: MPI_Comm_rank(A->comm,&rank);
226:
227: if (lu->StatPrint) { /* collect time for mat conversion */
228: MPI_Barrier(A->comm);
229: PetscGetTime(&time0);
230: }
232: if (lu->MatInputMode == GLOBAL) { /* global mat input */
233: if (size > 1) { /* convert mpi A to seq mat A */
234: ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);
235: MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);
236: ISDestroy(isrow);
237:
238: A_seq = *tseq;
239: PetscFree(tseq);
240: aa = (Mat_SeqAIJ*)A_seq->data;
241: } else {
242: aa = (Mat_SeqAIJ*)A->data;
243: }
245: /* Allocate storage, then convert Petsc NR matrix to SuperLU_DIST NC */
246: if (lu->flg == DIFFERENT_NONZERO_PATTERN) {/* first numeric factorization */
247: #if defined(PETSC_USE_COMPLEX)
248: zallocateA_dist(N, aa->nz, &lu->val, &lu->col, &lu->row);
249: #else
250: dallocateA_dist(N, aa->nz, &lu->val, &lu->col, &lu->row);
251: #endif
252: } else { /* successive numeric factorization, sparsity pattern is reused. */
253: Destroy_CompCol_Matrix_dist(&lu->A_sup);
254: Destroy_LU(N, &lu->grid, &lu->LUstruct);
255: lu->options.Fact = SamePattern;
256: }
257: #if defined(PETSC_USE_COMPLEX)
258: zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,&lu->val,&lu->col, &lu->row);
259: #else
260: dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,aa->j,aa->i,&lu->val, &lu->col, &lu->row);
261: #endif
263: /* Create compressed column matrix A_sup. */
264: #if defined(PETSC_USE_COMPLEX)
265: zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE);
266: #else
267: dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE);
268: #endif
269: } else { /* distributed mat input */
270: Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
271: aa=(Mat_SeqAIJ*)(mat->A)->data;
272: bb=(Mat_SeqAIJ*)(mat->B)->data;
273: ai=aa->i; aj=aa->j;
274: bi=bb->i; bj=bb->j;
275: #if defined(PETSC_USE_COMPLEX)
276: av=(doublecomplex*)aa->a;
277: bv=(doublecomplex*)bb->a;
278: #else
279: av=aa->a;
280: bv=bb->a;
281: #endif
282: rstart = mat->rstart;
283: nz = aa->nz + bb->nz;
284: garray = mat->garray;
285: rstart = mat->rstart;
287: if (lu->flg == DIFFERENT_NONZERO_PATTERN) {/* first numeric factorization */
288: #if defined(PETSC_USE_COMPLEX)
289: zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row);
290: #else
291: dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row);
292: #endif
293: } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
294: /* Destroy_CompRowLoc_Matrix_dist(&lu->A_sup); */ /* crash! */
295: Destroy_LU(N, &lu->grid, &lu->LUstruct);
296: lu->options.Fact = SamePattern;
297: }
298: nz = 0; irow = mat->rstart;
299: for ( i=0; i<m; i++ ) {
300: lu->row[i] = nz;
301: countA = ai[i+1] - ai[i];
302: countB = bi[i+1] - bi[i];
303: ajj = aj + ai[i]; /* ptr to the beginning of this row */
304: bjj = bj + bi[i];
306: /* B part, smaller col index */
307: colA_start = mat->rstart + ajj[0]; /* the smallest global col index of A */
308: jB = 0;
309: for (j=0; j<countB; j++){
310: jcol = garray[bjj[j]];
311: if (jcol > colA_start) {
312: jB = j;
313: break;
314: }
315: lu->col[nz] = jcol;
316: lu->val[nz++] = *bv++;
317: if (j==countB-1) jB = countB;
318: }
320: /* A part */
321: for (j=0; j<countA; j++){
322: lu->col[nz] = mat->rstart + ajj[j];
323: lu->val[nz++] = *av++;
324: }
326: /* B part, larger col index */
327: for (j=jB; j<countB; j++){
328: lu->col[nz] = garray[bjj[j]];
329: lu->val[nz++] = *bv++;
330: }
331: }
332: lu->row[m] = nz;
333: #if defined(PETSC_USE_COMPLEX)
334: zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,
335: lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE);
336: #else
337: dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,
338: lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE);
339: #endif
340: }
341: if (lu->StatPrint) {
342: PetscGetTime(&time);
343: time0 = time - time0;
344: }
346: /* Factor the matrix. */
347: PStatInit(&stat); /* Initialize the statistics variables. */
349: if (lu->MatInputMode == GLOBAL) { /* global mat input */
350: #if defined(PETSC_USE_COMPLEX)
351: pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,
352: &lu->grid, &lu->LUstruct, berr, &stat, &info);
353: #else
354: pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,
355: &lu->grid, &lu->LUstruct, berr, &stat, &info);
356: #endif
357: } else { /* distributed mat input */
358: #if defined(PETSC_USE_COMPLEX)
359: pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0, &lu->grid,
360: &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
361: if (info) SETERRQ1(1,"pzgssvx fails, info: %d\n",info);
362: #else
363: pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0, &lu->grid,
364: &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
365: if (info) SETERRQ1(1,"pdgssvx fails, info: %d\n",info);
366: #endif
367: }
369: if (lu->MatInputMode == GLOBAL && size > 1){
370: MatDestroy(A_seq);
371: }
373: if (lu->StatPrint) {
374: if (size > 1){
375: MPI_Reduce(&time0,&time_max,1,MPI_DOUBLE,MPI_MAX,0,A->comm);
376: MPI_Reduce(&time0,&time_min,1,MPI_DOUBLE,MPI_MIN,0,A->comm);
377: MPI_Reduce(&time0,&time,1,MPI_DOUBLE,MPI_SUM,0,A->comm);
378: time = time/size; /* average time */
379: if (!rank)
380: PetscPrintf(PETSC_COMM_SELF, " Mat conversion(PETSc->SuperLU_DIST) time (max/min/avg): \n \
381: %g / %g / %g\n",time_max,time_min,time);
382: } else {
383: PetscPrintf(PETSC_COMM_SELF, " Mat conversion(PETSc->SuperLU_DIST) time: \n \
384: %g\n",time0);
385: }
386:
387: PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */
388: }
389: PStatFree(&stat);
390: (*F)->assembled = PETSC_TRUE;
391: lu->flg = SAME_NONZERO_PATTERN;
392: return(0);
393: }
395: /* Note the Petsc r and c permutations are ignored */
398: int MatLUFactorSymbolic_SuperLU_DIST(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
399: {
400: Mat B;
401: Mat_SuperLU_DIST *lu;
402: int ierr,M=A->M,N=A->N,size,indx;
403: superlu_options_t options;
404: PetscTruth flg;
405: const char *ptype[] = {"MMD_AT_PLUS_A","NATURAL","MMD_ATA","COLAMD"};
406: const char *prtype[] = {"LargeDiag","NATURAL"};
409: /* Create the factorization matrix */
410: MatCreate(A->comm,A->m,A->n,M,N,&B);
411: MatSetType(B,A->type_name);
412: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
413: MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);
415: B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
416: B->ops->solve = MatSolve_SuperLU_DIST;
417: B->factor = FACTOR_LU;
419: lu = (Mat_SuperLU_DIST*)(B->spptr);
421: /* Set the input options */
422: set_default_options_dist(&options);
423: lu->MatInputMode = GLOBAL;
424: MPI_Comm_dup(A->comm,&(lu->comm_superlu));
426: MPI_Comm_size(A->comm,&size);
427: lu->nprow = size/2; /* Default process rows. */
428: if (lu->nprow == 0) lu->nprow = 1;
429: lu->npcol = size/lu->nprow; /* Default process columns. */
431: PetscOptionsBegin(A->comm,A->prefix,"SuperLU_Dist Options","Mat");
432:
433: PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,&lu->nprow,PETSC_NULL);
434: PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,&lu->npcol,PETSC_NULL);
435: if (size != lu->nprow * lu->npcol) SETERRQ(1,"Number of processes should be equal to nprow*npcol");
436:
437: PetscOptionsInt("-mat_superlu_dist_matinput","Matrix input mode (0: GLOBAL; 1: DISTRIBUTED)","None",lu->MatInputMode,&lu->MatInputMode,PETSC_NULL);
438: if(lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL;
440: PetscOptionsLogical("-mat_superlu_dist_equil","Equilibrate matrix","None",PETSC_TRUE,&flg,0);
441: if (!flg) {
442: options.Equil = NO;
443: }
445: PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",prtype,2,prtype[0],&indx,&flg);
446: if (flg) {
447: switch (indx) {
448: case 0:
449: options.RowPerm = LargeDiag;
450: break;
451: case 1:
452: options.RowPerm = NOROWPERM;
453: break;
454: }
455: }
457: PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",ptype,4,ptype[0],&indx,&flg);
458: if (flg) {
459: switch (indx) {
460: case 0:
461: options.ColPerm = MMD_AT_PLUS_A;
462: break;
463: case 1:
464: options.ColPerm = NATURAL;
465: break;
466: case 2:
467: options.ColPerm = MMD_ATA;
468: break;
469: case 3:
470: options.ColPerm = COLAMD;
471: break;
472: }
473: }
475: PetscOptionsLogical("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",PETSC_TRUE,&flg,0);
476: if (!flg) {
477: options.ReplaceTinyPivot = NO;
478: }
480: options.IterRefine = NOREFINE;
481: PetscOptionsLogical("-mat_superlu_dist_iterrefine","Use iterative refinement","None",PETSC_FALSE,&flg,0);
482: if (flg) {
483: options.IterRefine = DOUBLE;
484: }
486: if (PetscLogPrintInfo) {
487: lu->StatPrint = (int)PETSC_TRUE;
488: } else {
489: lu->StatPrint = (int)PETSC_FALSE;
490: }
491: PetscOptionsLogical("-mat_superlu_dist_statprint","Print factorization information","None",
492: (PetscTruth)lu->StatPrint,(PetscTruth*)&lu->StatPrint,0);
493: PetscOptionsEnd();
495: /* Initialize the SuperLU process grid. */
496: superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid);
498: /* Initialize ScalePermstruct and LUstruct. */
499: ScalePermstructInit(M, N, &lu->ScalePermstruct);
500: LUstructInit(M, N, &lu->LUstruct);
502: lu->options = options;
503: lu->flg = DIFFERENT_NONZERO_PATTERN;
504: lu->CleanUpSuperLU_Dist = PETSC_TRUE;
505: *F = B;
506: return(0);
507: }
511: int MatAssemblyEnd_SuperLU_DIST(Mat A,MatAssemblyType mode) {
512: int ierr;
513: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)(A->spptr);
516: (*lu->MatAssemblyEnd)(A,mode);
517: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
518: A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
519: return(0);
520: }
524: int MatFactorInfo_SuperLU_DIST(Mat A,PetscViewer viewer)
525: {
526: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)A->spptr;
527: superlu_options_t options;
528: int ierr;
529: char *colperm;
532: /* check if matrix is superlu_dist type */
533: if (A->ops->solve != MatSolve_SuperLU_DIST) return(0);
535: options = lu->options;
536: PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");
537: PetscViewerASCIIPrintf(viewer," Equilibrate matrix %s \n",(options.Equil != NO) ? "true": "false");
538: PetscViewerASCIIPrintf(viewer," Matrix input mode %d \n",lu->MatInputMode);
539: PetscViewerASCIIPrintf(viewer," Replace tiny pivots %s \n",(options.ReplaceTinyPivot != NO) ? "true": "false");
540: PetscViewerASCIIPrintf(viewer," Use iterative refinement %s \n",(options.IterRefine == DOUBLE) ? "true": "false");
541: PetscViewerASCIIPrintf(viewer," Processors in row %d col partition %d \n",lu->nprow,lu->npcol);
542: PetscViewerASCIIPrintf(viewer," Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL": "LargeDiag");
543: if (options.ColPerm == NATURAL) {
544: colperm = "NATURAL";
545: } else if (options.ColPerm == MMD_AT_PLUS_A) {
546: colperm = "MMD_AT_PLUS_A";
547: } else if (options.ColPerm == MMD_ATA) {
548: colperm = "MMD_ATA";
549: } else if (options.ColPerm == COLAMD) {
550: colperm = "COLAMD";
551: } else {
552: SETERRQ(1,"Unknown column permutation");
553: }
554: PetscViewerASCIIPrintf(viewer," Column permutation %s \n",colperm);
555: return(0);
556: }
560: int MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
561: {
562: int ierr;
563: PetscTruth isascii;
564: PetscViewerFormat format;
565: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)(A->spptr);
568: (*lu->MatView)(A,viewer);
570: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
571: if (isascii) {
572: PetscViewerGetFormat(viewer,&format);
573: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
574: MatFactorInfo_SuperLU_DIST(A,viewer);
575: }
576: }
577: return(0);
578: }
581: EXTERN_C_BEGIN
584: int MatConvert_Base_SuperLU_DIST(Mat A,const MatType type,Mat *newmat) {
585: /* This routine is only called to convert to MATSUPERLU_DIST */
586: /* from MATSEQAIJ if A has a single process communicator */
587: /* or MATMPIAIJ otherwise, so we will ignore 'MatType type'. */
588: int ierr, size;
589: MPI_Comm comm;
590: Mat B=*newmat;
591: Mat_SuperLU_DIST *lu;
594: if (B != A) {
595: MatDuplicate(A,MAT_COPY_VALUES,&B);
596: }
598: PetscObjectGetComm((PetscObject)A,&comm);
599: PetscNew(Mat_SuperLU_DIST,&lu);
601: lu->MatDuplicate = A->ops->duplicate;
602: lu->MatView = A->ops->view;
603: lu->MatAssemblyEnd = A->ops->assemblyend;
604: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
605: lu->MatDestroy = A->ops->destroy;
606: lu->CleanUpSuperLU_Dist = PETSC_FALSE;
608: B->spptr = (void*)lu;
609: B->ops->duplicate = MatDuplicate_SuperLU_DIST;
610: B->ops->view = MatView_SuperLU_DIST;
611: B->ops->assemblyend = MatAssemblyEnd_SuperLU_DIST;
612: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
613: B->ops->destroy = MatDestroy_SuperLU_DIST;
614: MPI_Comm_size(comm,&size);
615: if (size == 1) {
616: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_dist_C",
617: "MatConvert_Base_SuperLU_DIST",MatConvert_Base_SuperLU_DIST);
618: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_dist_seqaij_C",
619: "MatConvert_SuperLU_DIST_Base",MatConvert_SuperLU_DIST_Base);
620: } else {
621: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_superlu_dist_C",
622: "MatConvert_Base_SuperLU_DIST",MatConvert_Base_SuperLU_DIST);
623: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_dist_mpiaij_C",
624: "MatConvert_SuperLU_DIST_Base",MatConvert_SuperLU_DIST_Base);
625: }
626: PetscLogInfo(0,"Using SuperLU_DIST for SeqAIJ LU factorization and solves.");
627: PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU_DIST);
628: *newmat = B;
629: return(0);
630: }
631: EXTERN_C_END
635: int MatDuplicate_SuperLU_DIST(Mat A, MatDuplicateOption op, Mat *M) {
636: int ierr;
637: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST *)A->spptr;
640: (*lu->MatDuplicate)(A,op,M);
641: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU_DIST));
642: return(0);
643: }
645: /*MC
646: MATSUPERLU_DIST - MATSUPERLU_DIST = "superlu_dist" - A matrix type providing direct solvers (LU) for parallel matrices
647: via the external package SuperLU_DIST.
649: If SuperLU_DIST is installed (see the manual for
650: instructions on how to declare the existence of external packages),
651: a matrix type can be constructed which invokes SuperLU_DIST solvers.
652: After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU_DIST).
653: This matrix type is only supported for double precision real.
655: This matrix inherits from MATSEQAIJ when constructed with a single process communicator,
656: and from MATMPIAIJ otherwise. As a result, for single process communicators,
657: MatSeqAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
658: for communicators controlling multiple processes. It is recommended that you call both of
659: the above preallocation routines for simplicity. One can also call MatConvert for an inplace
660: conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
661: without data copy.
663: Options Database Keys:
664: + -mat_type superlu_dist - sets the matrix type to "superlu_dist" during a call to MatSetFromOptions()
665: . -mat_superlu_dist_r <n> - number of rows in processor partition
666: . -mat_superlu_dist_c <n> - number of columns in processor partition
667: . -mat_superlu_dist_matinput <0,1> - matrix input mode; 0=global, 1=distributed
668: . -mat_superlu_dist_equil - equilibrate the matrix
669: . -mat_superlu_dist_rowperm <LargeDiag,NATURAL> - row permutation
670: . -mat_superlu_dist_colperm <MMD_AT_PLUS_A,MMD_ATA,COLAMD,NATURAL> - column permutation
671: . -mat_superlu_dist_replacetinypivot - replace tiny pivots
672: . -mat_superlu_dist_iterrefine - use iterative refinement
673: - -mat_superlu_dist_statprint - print factorization information
675: Level: beginner
677: .seealso: PCLU
678: M*/
680: EXTERN_C_BEGIN
683: int MatCreate_SuperLU_DIST(Mat A) {
684: int ierr,size;
685: Mat A_diag;
688: /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
689: /* and SuperLU_DIST types */
690: PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU_DIST);
691: MPI_Comm_size(A->comm,&size);
692: if (size == 1) {
693: MatSetType(A,MATSEQAIJ);
694: } else {
695: MatSetType(A,MATMPIAIJ);
696: A_diag = ((Mat_MPIAIJ *)A->data)->A;
697: MatConvert_Base_SuperLU_DIST(A_diag,MATSUPERLU_DIST,&A_diag);
698: }
699: MatConvert_Base_SuperLU_DIST(A,MATSUPERLU_DIST,&A);
700: return(0);
701: }
702: EXTERN_C_END