Actual source code: superlu_dist.c

petsc-master 2016-12-01
Report Typos and Errors

  2: /*
  3:         Provides an interface to the SuperLU_DIST 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 around weird problem with SuperLU on cray */
  9: #include <stdlib.h>
 10: #endif

 12: #if defined(PETSC_USE_64BIT_INDICES)
 13: /* ugly SuperLU_Dist variable telling it to use long long int */
 14: #define _LONGINT
 15: #endif

 17: EXTERN_C_BEGIN
 18: #if defined(PETSC_USE_COMPLEX)
 19: #include <superlu_zdefs.h>
 20: #else
 21: #include <superlu_ddefs.h>
 22: #endif
 23: EXTERN_C_END

 25: /*
 26:     GLOBAL - The sparse matrix and right hand side are all stored initially on process 0. Should be called centralized
 27:     DISTRIBUTED - The sparse matrix and right hand size are initially stored across the entire MPI communicator.
 28: */
 29: typedef enum {GLOBAL,DISTRIBUTED} SuperLU_MatInputMode;
 30: const char *SuperLU_MatInputModes[] = {"GLOBAL","DISTRIBUTED","SuperLU_MatInputMode","PETSC_",0};

 32: typedef struct {
 33:   int_t                  nprow,npcol,*row,*col;
 34:   gridinfo_t             grid;
 35:   superlu_dist_options_t options;
 36:   SuperMatrix            A_sup;
 37:   ScalePermstruct_t      ScalePermstruct;
 38:   LUstruct_t             LUstruct;
 39:   int                    StatPrint;
 40:   SuperLU_MatInputMode   MatInputMode;
 41:   SOLVEstruct_t          SOLVEstruct;
 42:   fact_t                 FactPattern;
 43:   MPI_Comm               comm_superlu;
 44: #if defined(PETSC_USE_COMPLEX)
 45:   doublecomplex          *val;
 46: #else
 47:   double                 *val;
 48: #endif
 49:   PetscBool              matsolve_iscalled,matmatsolve_iscalled;
 50:   PetscBool              CleanUpSuperLU_Dist;  /* Flag to clean up (non-global) SuperLU objects during Destroy */
 51: } Mat_SuperLU_DIST;


 56: PetscErrorCode MatSuperluDistGetDiagU_SuperLU_DIST(Mat F,PetscScalar *diagU)
 57: {
 58:   Mat_SuperLU_DIST  *lu= (Mat_SuperLU_DIST*)F->data;

 61: #if defined(PETSC_USE_COMPLEX)
 62:   PetscStackCall("SuperLU_DIST:pzGetDiagU",pzGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,(doublecomplex*)diagU));
 63: #else
 64:   PetscStackCall("SuperLU_DIST:pdGetDiagU",pdGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,diagU));
 65: #endif
 66:   return(0);
 67: }

 71: PetscErrorCode MatSuperluDistGetDiagU(Mat F,PetscScalar *diagU)
 72: {
 73:   PetscErrorCode    ierr;

 77:   PetscTryMethod(F,"MatSuperluDistGetDiagU_C",(Mat,PetscScalar*),(F,diagU));
 78:   return(0);
 79: }

 83: static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
 84: {
 85:   PetscErrorCode   ierr;
 86:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;

 89:   if (lu->CleanUpSuperLU_Dist) {
 90:     /* Deallocate SuperLU_DIST storage */
 91:     if (lu->MatInputMode == GLOBAL) {
 92:       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
 93:     } else {
 94:       PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
 95:       if (lu->options.SolveInitialized) {
 96: #if defined(PETSC_USE_COMPLEX)
 97:         PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
 98: #else
 99:         PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
100: #endif
101:       }
102:     }
103:     PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
104:     PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
105:     PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));

107:     /* Release the SuperLU_DIST process grid. */
108:     PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
109:     MPI_Comm_free(&(lu->comm_superlu));
110:   }
111:   PetscFree(A->data);
112:   /* clear composed functions */
113:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
114:   PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);

116:   return(0);
117: }

121: static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
122: {
123:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
124:   PetscErrorCode   ierr;
125:   PetscMPIInt      size;
126:   PetscInt         m=A->rmap->n,M=A->rmap->N,N=A->cmap->N;
127:   SuperLUStat_t    stat;
128:   double           berr[1];
129:   PetscScalar      *bptr;
130:   PetscInt         nrhs=1;
131:   Vec              x_seq;
132:   IS               iden;
133:   VecScatter       scat;
134:   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
135:   static PetscBool cite = PETSC_FALSE;

138:   PetscCitationsRegister("@article{lidemmel03,\n  author = {Xiaoye S. Li and James W. Demmel},\n  title = {{SuperLU_DIST}: A Scalable Distributed-Memory Sparse Direct\n           Solver for Unsymmetric Linear Systems},\n  journal = {ACM Trans. Mathematical Software},\n  volume = {29},\n  number = {2},\n  pages = {110-140},\n  year = 2003\n}\n",&cite);
139:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
140:   if (size > 1 && lu->MatInputMode == GLOBAL) {
141:     /* global mat input, convert b to x_seq */
142:     VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);
143:     ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);
144:     VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);
145:     ISDestroy(&iden);

147:     VecScatterBegin(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);
148:     VecScatterEnd(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);
149:     VecGetArray(x_seq,&bptr);
150:   } else { /* size==1 || distributed mat input */
151:     if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
152:       /* see comments in MatMatSolve() */
153: #if defined(PETSC_USE_COMPLEX)
154:       PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
155: #else
156:       PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
157: #endif
158:       lu->options.SolveInitialized = NO;
159:     }
160:     VecCopy(b_mpi,x);
161:     VecGetArray(x,&bptr);
162:   }

164:   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");

166:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
167:   if (lu->MatInputMode == GLOBAL) {
168: #if defined(PETSC_USE_COMPLEX)
169:     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
170: #else
171:     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
172: #endif
173:   } else { /* distributed mat input */
174: #if defined(PETSC_USE_COMPLEX)
175:     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
176: #else
177:     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
178: #endif
179:   }
180:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);

182:   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid);      /* Print the statistics. */
183:   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));

185:   if (size > 1 && lu->MatInputMode == GLOBAL) {
186:     /* convert seq x to mpi x */
187:     VecRestoreArray(x_seq,&bptr);
188:     VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);
189:     VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);
190:     VecScatterDestroy(&scat);
191:     VecDestroy(&x_seq);
192:   } else {
193:     VecRestoreArray(x,&bptr);

195:     lu->matsolve_iscalled    = PETSC_TRUE;
196:     lu->matmatsolve_iscalled = PETSC_FALSE;
197:   }
198:   return(0);
199: }

203: static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
204: {
205:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
206:   PetscErrorCode   ierr;
207:   PetscMPIInt      size;
208:   PetscInt         M=A->rmap->N,m=A->rmap->n,nrhs;
209:   SuperLUStat_t    stat;
210:   double           berr[1];
211:   PetscScalar      *bptr;
212:   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
213:   PetscBool        flg;

216:   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");
217:   PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
218:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
219:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
220:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");

222:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
223:   if (size > 1 && lu->MatInputMode == GLOBAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatInputMode=GLOBAL for nproc %d>1 is not supported",size);
224:   /* size==1 or distributed mat input */
225:   if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
226:     /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
227:        thus destroy it and create a new SOLVEstruct.
228:        Otherwise it may result in memory corruption or incorrect solution
229:        See src/mat/examples/tests/ex125.c */
230: #if defined(PETSC_USE_COMPLEX)
231:     PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
232: #else
233:     PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
234: #endif
235:     lu->options.SolveInitialized = NO;
236:   }
237:   MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);

239:   MatGetSize(B_mpi,NULL,&nrhs);

241:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
242:   MatDenseGetArray(X,&bptr);
243:   if (lu->MatInputMode == GLOBAL) { /* size == 1 */
244: #if defined(PETSC_USE_COMPLEX)
245:     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, M, nrhs,&lu->grid, &lu->LUstruct, berr, &stat, &info));
246: #else
247:     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, M, nrhs, &lu->grid, &lu->LUstruct, berr, &stat, &info));
248: #endif
249:   } else { /* distributed mat input */
250: #if defined(PETSC_USE_COMPLEX)
251:     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid, &lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
252: #else
253:     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
254: #endif
255:   }
256:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
257:   MatDenseRestoreArray(X,&bptr);

259:   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */
260:   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
261:   lu->matsolve_iscalled    = PETSC_FALSE;
262:   lu->matmatsolve_iscalled = PETSC_TRUE;
263:   return(0);
264: }


269: static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
270: {
271:   Mat              *tseq,A_seq = NULL;
272:   Mat_SeqAIJ       *aa,*bb;
273:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)(F)->data;
274:   PetscErrorCode   ierr;
275:   PetscInt         M=A->rmap->N,N=A->cmap->N,i,*ai,*aj,*bi,*bj,nz,rstart,*garray,
276:                    m=A->rmap->n, colA_start,j,jcol,jB,countA,countB,*bjj,*ajj;
277:   int              sinfo;   /* SuperLU_Dist info flag is always an int even with long long indices */
278:   PetscMPIInt      size;
279:   SuperLUStat_t    stat;
280:   double           *berr=0;
281:   IS               isrow;
282: #if defined(PETSC_USE_COMPLEX)
283:   doublecomplex    *av, *bv;
284: #else
285:   double           *av, *bv;
286: #endif

289:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);

291:   if (lu->MatInputMode == GLOBAL) { /* global mat input */
292:     if (size > 1) { /* convert mpi A to seq mat A */
293:       ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);
294:       MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);
295:       ISDestroy(&isrow);

297:       A_seq = *tseq;
298:       PetscFree(tseq);
299:       aa    = (Mat_SeqAIJ*)A_seq->data;
300:     } else {
301:       PetscBool flg;
302:       PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);
303:       if (flg) {
304:         Mat_MPIAIJ *At = (Mat_MPIAIJ*)A->data;
305:         A = At->A;
306:       }
307:       aa =  (Mat_SeqAIJ*)A->data;
308:     }

310:     /* Convert Petsc NR matrix to SuperLU_DIST NC.
311:        Note: memories of lu->val, col and row are allocated by CompRow_to_CompCol_dist()! */
312:     if (lu->options.Fact != DOFACT) {/* successive numeric factorization, sparsity pattern is reused. */
313:       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
314:       if (lu->FactPattern == SamePattern_SameRowPerm) {
315:         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
316:       } else { /* lu->FactPattern == SamePattern */
317:         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct));
318:         lu->options.Fact = SamePattern;
319:       }
320:     }
321: #if defined(PETSC_USE_COMPLEX)
322:     PetscStackCall("SuperLU_DIST:zCompRow_to_CompCol_dist",zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val,&lu->col, &lu->row));
323: #else
324:     PetscStackCall("SuperLU_DIST:dCompRow_to_CompCol_dist",dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val, &lu->col, &lu->row));
325: #endif

327:     /* Create compressed column matrix A_sup. */
328: #if defined(PETSC_USE_COMPLEX)
329:     PetscStackCall("SuperLU_DIST:zCreate_CompCol_Matrix_dist",zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE));
330: #else
331:     PetscStackCall("SuperLU_DIST:dCreate_CompCol_Matrix_dist",dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE));
332: #endif
333:   } else { /* distributed mat input */
334:     Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
335:     aa=(Mat_SeqAIJ*)(mat->A)->data;
336:     bb=(Mat_SeqAIJ*)(mat->B)->data;
337:     ai=aa->i; aj=aa->j;
338:     bi=bb->i; bj=bb->j;
339: #if defined(PETSC_USE_COMPLEX)
340:     av=(doublecomplex*)aa->a;
341:     bv=(doublecomplex*)bb->a;
342: #else
343:     av=aa->a;
344:     bv=bb->a;
345: #endif
346:     rstart = A->rmap->rstart;
347:     nz     = aa->nz + bb->nz;
348:     garray = mat->garray;

350:     if (lu->options.Fact == DOFACT) { /* first numeric factorization */
351: #if defined(PETSC_USE_COMPLEX)
352:       PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
353: #else
354:       PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
355: #endif
356:     } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
357:       if (lu->FactPattern == SamePattern_SameRowPerm) {
358:         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
359:       } else if (lu->FactPattern == SamePattern) {
360:         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); /* Deallocate storage associated with the L and U matrices. */
361:         lu->options.Fact = SamePattern;
362:       } else {
363:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm");
364:       }
365:     }
366:     nz = 0;
367:     for (i=0; i<m; i++) {
368:       lu->row[i] = nz;
369:       countA     = ai[i+1] - ai[i];
370:       countB     = bi[i+1] - bi[i];
371:       if (aj) {
372:         ajj = aj + ai[i]; /* ptr to the beginning of this row */
373:       } else {
374:         ajj = NULL;
375:       }
376:       bjj = bj + bi[i];

378:       /* B part, smaller col index */
379:       if (aj) {
380:         colA_start = rstart + ajj[0]; /* the smallest global col index of A */
381:       } else { /* superlu_dist does not require matrix has diagonal entries, thus aj=NULL would work */
382:         colA_start = rstart;
383:       }
384:       jB         = 0;
385:       for (j=0; j<countB; j++) {
386:         jcol = garray[bjj[j]];
387:         if (jcol > colA_start) {
388:           jB = j;
389:           break;
390:         }
391:         lu->col[nz]   = jcol;
392:         lu->val[nz++] = *bv++;
393:         if (j==countB-1) jB = countB;
394:       }

396:       /* A part */
397:       for (j=0; j<countA; j++) {
398:         lu->col[nz]   = rstart + ajj[j];
399:         lu->val[nz++] = *av++;
400:       }

402:       /* B part, larger col index */
403:       for (j=jB; j<countB; j++) {
404:         lu->col[nz]   = garray[bjj[j]];
405:         lu->val[nz++] = *bv++;
406:       }
407:     }
408:     lu->row[m] = nz;

410:     if (lu->options.Fact == DOFACT) {
411: #if defined(PETSC_USE_COMPLEX)
412:       PetscStackCall("SuperLU_DIST:zCreate_CompRowLoc_Matrix_dist",zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE));
413: #else
414:       PetscStackCall("SuperLU_DIST:dCreate_CompRowLoc_Matrix_dist",dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE));
415: #endif
416:     }
417:   }

419:   /* Factor the matrix. */
420:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));   /* Initialize the statistics variables. */
421:   if (lu->MatInputMode == GLOBAL) { /* global mat input */
422: #if defined(PETSC_USE_COMPLEX)
423:     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
424: #else
425:     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
426: #endif
427:   } else { /* distributed mat input */
428: #if defined(PETSC_USE_COMPLEX)
429:     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
430: #else
431:     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
432: #endif
433:   }
434: 
435:   if (sinfo > 0) {
436:     if (A->erroriffailure) {
437:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
438:     } else {
439:       if (sinfo <= lu->A_sup.ncol) {
440:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
441:         PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);
442:       } else if (sinfo > lu->A_sup.ncol) {
443:         /* 
444:          number of bytes allocated when memory allocation
445:          failure occurred, plus A->ncol.
446:          */
447:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
448:         PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);
449:       }
450:     }
451:   } else if (sinfo < 0) {
452:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo);
453:   }

455:   if (lu->MatInputMode == GLOBAL && size > 1) {
456:     MatDestroy(&A_seq);
457:   }

459:   if (lu->options.PrintStat) {
460:     PStatPrint(&lu->options, &stat, &lu->grid);  /* Print the statistics. */
461:   }
462:   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
463:   (F)->assembled    = PETSC_TRUE;
464:   (F)->preallocated = PETSC_TRUE;
465:   lu->options.Fact  = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
466:   return(0);
467: }

469: /* Note the Petsc r and c permutations are ignored */
472: static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
473: {
474:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
475:   PetscInt         M   = A->rmap->N,N=A->cmap->N;

478:   /* Initialize the SuperLU process grid. */
479:   PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid));

481:   /* Initialize ScalePermstruct and LUstruct. */
482:   PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
483:   PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
484:   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
485:   F->ops->solve           = MatSolve_SuperLU_DIST;
486:   F->ops->matsolve        = MatMatSolve_SuperLU_DIST;
487:   lu->CleanUpSuperLU_Dist = PETSC_TRUE;
488:   return(0);
489: }

493: static PetscErrorCode MatFactorGetSolverPackage_aij_superlu_dist(Mat A,const MatSolverPackage *type)
494: {
496:   *type = MATSOLVERSUPERLU_DIST;
497:   return(0);
498: }

502: static PetscErrorCode MatFactorInfo_SuperLU_DIST(Mat A,PetscViewer viewer)
503: {
504:   Mat_SuperLU_DIST       *lu=(Mat_SuperLU_DIST*)A->data;
505:   superlu_dist_options_t options;
506:   PetscErrorCode         ierr;

509:   /* check if matrix is superlu_dist type */
510:   if (A->ops->solve != MatSolve_SuperLU_DIST) return(0);

512:   options = lu->options;
513:   PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");
514:   PetscViewerASCIIPrintf(viewer,"  Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);
515:   PetscViewerASCIIPrintf(viewer,"  Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);
516:   PetscViewerASCIIPrintf(viewer,"  Matrix input mode %d \n",lu->MatInputMode);
517:   PetscViewerASCIIPrintf(viewer,"  Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);
518:   PetscViewerASCIIPrintf(viewer,"  Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);
519:   PetscViewerASCIIPrintf(viewer,"  Processors in row %d col partition %d \n",lu->nprow,lu->npcol);
520:   PetscViewerASCIIPrintf(viewer,"  Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL" : "LargeDiag");

522:   switch (options.ColPerm) {
523:   case NATURAL:
524:     PetscViewerASCIIPrintf(viewer,"  Column permutation NATURAL\n");
525:     break;
526:   case MMD_AT_PLUS_A:
527:     PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_AT_PLUS_A\n");
528:     break;
529:   case MMD_ATA:
530:     PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_ATA\n");
531:     break;
532:   case METIS_AT_PLUS_A:
533:     PetscViewerASCIIPrintf(viewer,"  Column permutation METIS_AT_PLUS_A\n");
534:     break;
535:   case PARMETIS:
536:     PetscViewerASCIIPrintf(viewer,"  Column permutation PARMETIS\n");
537:     break;
538:   default:
539:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
540:   }

542:   PetscViewerASCIIPrintf(viewer,"  Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);

544:   if (lu->FactPattern == SamePattern) {
545:     PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern\n");
546:   } else {
547:     PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern_SameRowPerm\n");
548:   }
549:   return(0);
550: }

554: static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
555: {
556:   PetscErrorCode    ierr;
557:   PetscBool         iascii;
558:   PetscViewerFormat format;

561:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
562:   if (iascii) {
563:     PetscViewerGetFormat(viewer,&format);
564:     if (format == PETSC_VIEWER_ASCII_INFO) {
565:       MatFactorInfo_SuperLU_DIST(A,viewer);
566:     }
567:   }
568:   return(0);
569: }

573: static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
574: {
575:   Mat                    B;
576:   Mat_SuperLU_DIST       *lu;
577:   PetscErrorCode         ierr;
578:   PetscInt               M=A->rmap->N,N=A->cmap->N,indx;
579:   PetscMPIInt            size;
580:   superlu_dist_options_t options;
581:   PetscBool              flg;
582:   const char             *colperm[]     = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
583:   const char             *rowperm[]     = {"LargeDiag","NATURAL"};
584:   const char             *factPattern[] = {"SamePattern","SamePattern_SameRowPerm"};
585:   PetscBool              set;

588:   /* Create the factorization matrix */
589:   MatCreate(PetscObjectComm((PetscObject)A),&B);
590:   MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);
591:   PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);
592:   MatSetUp(B);
593:   B->ops->getinfo          = MatGetInfo_External;
594:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
595:   B->ops->view             = MatView_SuperLU_DIST;
596:   B->ops->destroy          = MatDestroy_SuperLU_DIST;

598:   B->factortype = MAT_FACTOR_LU;

600:   /* set solvertype */
601:   PetscFree(B->solvertype);
602:   PetscStrallocpy(MATSOLVERSUPERLU_DIST,&B->solvertype);

604:   PetscNewLog(B,&lu);
605:   B->data = lu;

607:   /* Set the default input options:
608:      options.Fact              = DOFACT;
609:      options.Equil             = YES;
610:      options.ParSymbFact       = NO;
611:      options.ColPerm           = METIS_AT_PLUS_A;
612:      options.RowPerm           = LargeDiag;
613:      options.ReplaceTinyPivot  = YES;
614:      options.IterRefine        = DOUBLE;
615:      options.Trans             = NOTRANS;
616:      options.SolveInitialized  = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
617:      options.RefineInitialized = NO;
618:      options.PrintStat         = YES;
619:   */
620:   set_default_options_dist(&options);

622:   MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(lu->comm_superlu));
623:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
624:   /* Default num of process columns and rows */
625:   lu->npcol = (int_t) (0.5 + PetscSqrtReal((PetscReal)size));
626:   if (!lu->npcol) lu->npcol = 1;
627:   while (lu->npcol > 0) {
628:     lu->nprow = (int_t) (size/lu->npcol);
629:     if (size == lu->nprow * lu->npcol) break;
630:     lu->npcol--;
631:   }

633:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");
634:   PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);
635:   PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);
636:   if (size != lu->nprow * lu->npcol) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of processes %d must equal to nprow %d * npcol %d",size,lu->nprow,lu->npcol);

638:   lu->MatInputMode = DISTRIBUTED;

640:   PetscOptionsEnum("-mat_superlu_dist_matinput","Matrix input mode (global or distributed)","None",SuperLU_MatInputModes,(PetscEnum)lu->MatInputMode,(PetscEnum*)&lu->MatInputMode,NULL);
641:   if (lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL;

643:   PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
644:   if (set && !flg) options.Equil = NO;

646:   PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,2,rowperm[0],&indx,&flg);
647:   if (flg) {
648:     switch (indx) {
649:     case 0:
650:       options.RowPerm = LargeDiag;
651:       break;
652:     case 1:
653:       options.RowPerm = NOROWPERM;
654:       break;
655:     }
656:   }

658:   PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);
659:   if (flg) {
660:     switch (indx) {
661:     case 0:
662:       options.ColPerm = NATURAL;
663:       break;
664:     case 1:
665:       options.ColPerm = MMD_AT_PLUS_A;
666:       break;
667:     case 2:
668:       options.ColPerm = MMD_ATA;
669:       break;
670:     case 3:
671:       options.ColPerm = METIS_AT_PLUS_A;
672:       break;
673:     case 4:
674:       options.ColPerm = PARMETIS;   /* only works for np>1 */
675:       break;
676:     default:
677:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
678:     }
679:   }

681:   options.ReplaceTinyPivot = NO;
682:   PetscOptionsBool("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",options.ReplaceTinyPivot ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
683:   if (set && flg) options.ReplaceTinyPivot = YES;

685:   options.ParSymbFact = NO;
686:   PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);
687:   if (set && flg && size>1) {
688:     if (lu->MatInputMode == GLOBAL) {
689: #if defined(PETSC_USE_INFO)
690:       PetscInfo(A,"Warning: '-mat_superlu_dist_parsymbfact' is ignored because MatInputMode=GLOBAL\n");
691: #endif
692:     } else {
693: #if defined(PETSC_HAVE_PARMETIS)
694:       options.ParSymbFact = YES;
695:       options.ColPerm     = PARMETIS;   /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */
696: #else
697:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS");
698: #endif
699:     }
700:   }

702:   lu->FactPattern = SamePattern_SameRowPerm;
703:   PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,2,factPattern[0],&indx,&flg);
704:   if (flg) {
705:     switch (indx) {
706:     case 0:
707:       lu->FactPattern = SamePattern;
708:       break;
709:     case 1:
710:       lu->FactPattern = SamePattern_SameRowPerm;
711:       break;
712:     }
713:   }

715:   options.IterRefine = NOREFINE;
716:   PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);
717:   if (set) {
718:     if (flg) options.IterRefine = SLU_DOUBLE;
719:     else options.IterRefine = NOREFINE;
720:   }

722:   if (PetscLogPrintInfo) options.PrintStat = YES;
723:   else options.PrintStat = NO;
724:   PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);
725:   PetscOptionsEnd();

727:   lu->options              = options;
728:   lu->options.Fact         = DOFACT;
729:   lu->matsolve_iscalled    = PETSC_FALSE;
730:   lu->matmatsolve_iscalled = PETSC_FALSE;

732:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_superlu_dist);
733:   PetscObjectComposeFunction((PetscObject)B,"MatSuperluDistGetDiagU_C",MatSuperluDistGetDiagU_SuperLU_DIST);

735:   *F = B;
736:   return(0);
737: }

741: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU_DIST(void)
742: {
745:   MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
746:   MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
747:   return(0);
748: }

750: /*MC
751:   MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization

753:   Use ./configure --download-superlu_dist --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with SuperLU_DIST

755:   Use -pc_type lu -pc_factor_mat_solver_package superlu_dist to us this direct solver

757:    Works with AIJ matrices

759:   Options Database Keys:
760: + -mat_superlu_dist_r <n> - number of rows in processor partition
761: . -mat_superlu_dist_c <n> - number of columns in processor partition
762: . -mat_superlu_dist_matinput <0,1> - matrix input mode; 0=global, 1=distributed
763: . -mat_superlu_dist_equil - equilibrate the matrix
764: . -mat_superlu_dist_rowperm <LargeDiag,NATURAL> - row permutation
765: . -mat_superlu_dist_colperm <MMD_AT_PLUS_A,MMD_ATA,NATURAL> - column permutation
766: . -mat_superlu_dist_replacetinypivot - replace tiny pivots
767: . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm
768: . -mat_superlu_dist_iterrefine - use iterative refinement
769: - -mat_superlu_dist_statprint - print factorization information

771:    Level: beginner

773: .seealso: PCLU

775: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage

777: M*/