Actual source code: superlu_dist.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: /*
  2:         Provides an interface to the SuperLU_DIST sparse solver
  3: */

  5: #include <../src/mat/impls/aij/seq/aij.h>
  6: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  7: #include <petscpkg_version.h>

  9: EXTERN_C_BEGIN
 10: #if defined(PETSC_USE_COMPLEX)
 11: #include <superlu_zdefs.h>
 12: #if PETSC_PKG_SUPERLU_DIST_VERSION_GE(6,3,0)
 13: #define LUstructInit zLUstructInit
 14: #define ScalePermstructInit zScalePermstructInit
 15: #define ScalePermstructFree zScalePermstructFree
 16: #define LUstructFree zLUstructFree
 17: #define Destroy_LU zDestroy_LU
 18: #define ScalePermstruct_t zScalePermstruct_t
 19: #define LUstruct_t zLUstruct_t
 20: #define SOLVEstruct_t zSOLVEstruct_t
 21: #endif
 22: #else
 23: #include <superlu_ddefs.h>
 24: #if PETSC_PKG_SUPERLU_DIST_VERSION_GE(6,3,0)
 25: #define LUstructInit dLUstructInit
 26: #define ScalePermstructInit dScalePermstructInit
 27: #define ScalePermstructFree dScalePermstructFree
 28: #define LUstructFree dLUstructFree
 29: #define Destroy_LU dDestroy_LU
 30: #define ScalePermstruct_t dScalePermstruct_t
 31: #define LUstruct_t dLUstruct_t
 32: #define SOLVEstruct_t dSOLVEstruct_t
 33: #endif
 34: #endif
 35: EXTERN_C_END

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


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

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

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

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

 81: /*  This allows reusing the Superlu_DIST communicator and grid when only a single SuperLU_DIST matrix is used at a time */
 82: typedef struct {
 83:   MPI_Comm   comm;
 84:   PetscBool  busy;
 85:   gridinfo_t grid;
 86: } PetscSuperLU_DIST;
 87: static PetscMPIInt Petsc_Superlu_dist_keyval = MPI_KEYVAL_INVALID;

 89: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Superlu_dist_keyval_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
 90: {
 91:   PetscErrorCode    ierr;
 92:   PetscSuperLU_DIST *context = (PetscSuperLU_DIST *) attr_val;

 95:   if (keyval != Petsc_Superlu_dist_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
 96:   PetscInfo(NULL,"Removing Petsc_Superlu_dist_keyval attribute from communicator that is being freed\n");
 97:   PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&context->grid));
 98:   MPI_Comm_free(&context->comm);
 99:   PetscFree(context);
100:   PetscFunctionReturn(MPI_SUCCESS);
101: }

103: /*
104:    Performs MPI_Comm_free_keyval() on Petsc_Superlu_dist_keyval but keeps the global variable for
105:    users who do not destroy all PETSc objects before PetscFinalize().

107:    The value Petsc_Superlu_dist_keyval is retained so that Petsc_Superlu_dist_keyval_Delete_Fn()
108:    can still check that the keyval associated with the MPI communicator is correct when the MPI
109:    communicator is destroyed.

111:    This is called in PetscFinalize()
112: */
113: static PetscErrorCode Petsc_Superlu_dist_keyval_free(void)
114: {
116:   PetscMPIInt    Petsc_Superlu_dist_keyval_temp = Petsc_Superlu_dist_keyval;

119:   PetscInfo(NULL,"Freeing Petsc_Superlu_dist_keyval\n");
120:   MPI_Comm_free_keyval(&Petsc_Superlu_dist_keyval_temp);
121:   return(0);
122: }

124: static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
125: {
126:   PetscErrorCode   ierr;
127:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;

130:   if (lu->CleanUpSuperLU_Dist) {
131:     /* Deallocate SuperLU_DIST storage */
132:     PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
133:     if (lu->options.SolveInitialized) {
134: #if defined(PETSC_USE_COMPLEX)
135:       PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
136: #else
137:       PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
138: #endif
139:     }
140:     PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
141:     PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
142:     PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));

144:     /* Release the SuperLU_DIST process grid. Only if the matrix has its own copy, this is it is not in the communicator context */
145:     if (lu->comm_superlu) {
146:       PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
147:       MPI_Comm_free(&(lu->comm_superlu));
148:     } else {
149:       PetscSuperLU_DIST *context;
150:       MPI_Comm          comm;
151:       PetscMPIInt       flg;

153:       PetscObjectGetComm((PetscObject)A,&comm);
154:       MPI_Comm_get_attr(comm,Petsc_Superlu_dist_keyval,&context,&flg);
155:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Communicator does not have expected Petsc_Superlu_dist_keyval attribute");
156:       context->busy = PETSC_FALSE;
157:     }
158:   }
159:   PetscFree(A->data);
160:   /* clear composed functions */
161:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
162:   PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);

164:   return(0);
165: }

167: static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
168: {
169:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
170:   PetscErrorCode   ierr;
171:   PetscInt         m=A->rmap->n;
172:   SuperLUStat_t    stat;
173:   double           berr[1];
174:   PetscScalar      *bptr=NULL;
175:   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
176:   static PetscBool cite = PETSC_FALSE;

179:   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
180:   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);

182:   if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
183:     /* see comments in MatMatSolve() */
184: #if defined(PETSC_USE_COMPLEX)
185:     PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
186: #else
187:     PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
188: #endif
189:     lu->options.SolveInitialized = NO;
190:   }
191:   VecCopy(b_mpi,x);
192:   VecGetArray(x,&bptr);

194:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
195: #if defined(PETSC_USE_COMPLEX)
196:   PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,1,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
197: #else
198:   PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,1,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
199: #endif
200:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);

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

205:   VecRestoreArray(x,&bptr);
206:   lu->matsolve_iscalled    = PETSC_TRUE;
207:   lu->matmatsolve_iscalled = PETSC_FALSE;
208:   return(0);
209: }

211: static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
212: {
213:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
214:   PetscErrorCode   ierr;
215:   PetscInt         m=A->rmap->n,nrhs;
216:   SuperLUStat_t    stat;
217:   double           berr[1];
218:   PetscScalar      *bptr;
219:   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
220:   PetscBool        flg;

223:   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
224:   PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
225:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
226:   if (X != B_mpi) {
227:     PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
228:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
229:   }

231:   if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
232:     /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
233:        thus destroy it and create a new SOLVEstruct.
234:        Otherwise it may result in memory corruption or incorrect solution
235:        See src/mat/tests/ex125.c */
236: #if defined(PETSC_USE_COMPLEX)
237:     PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
238: #else
239:     PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
240: #endif
241:     lu->options.SolveInitialized = NO;
242:   }
243:   if (X != B_mpi) {
244:     MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);
245:   }

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

249:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
250:   MatDenseGetArray(X,&bptr);

252: #if defined(PETSC_USE_COMPLEX)
253:   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));
254: #else
255:   PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
256: #endif

258:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
259:   MatDenseRestoreArray(X,&bptr);

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

268: /*
269:   input:
270:    F:        numeric Cholesky factor
271:   output:
272:    nneg:     total number of negative pivots
273:    nzero:    total number of zero pivots
274:    npos:     (global dimension of F) - nneg - nzero
275: */
276: static PetscErrorCode MatGetInertia_SuperLU_DIST(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
277: {
278:   PetscErrorCode   ierr;
279:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
280:   PetscScalar      *diagU=NULL;
281:   PetscInt         M,i,neg=0,zero=0,pos=0;
282:   PetscReal        r;

285:   if (!F->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix factor F is not assembled");
286:   if (lu->options.RowPerm != NOROWPERM) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must set NOROWPERM");
287:   MatGetSize(F,&M,NULL);
288:   PetscMalloc1(M,&diagU);
289:   MatSuperluDistGetDiagU(F,diagU);
290:   for (i=0; i<M; i++) {
291: #if defined(PETSC_USE_COMPLEX)
292:     r = PetscImaginaryPart(diagU[i])/10.0;
293:     if (r< -PETSC_MACHINE_EPSILON || r>PETSC_MACHINE_EPSILON) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"diagU[%d]=%g + i %g is non-real",i,PetscRealPart(diagU[i]),r*10.0);
294:     r = PetscRealPart(diagU[i]);
295: #else
296:     r = diagU[i];
297: #endif
298:     if (r > 0) {
299:       pos++;
300:     } else if (r < 0) {
301:       neg++;
302:     } else zero++;
303:   }

305:   PetscFree(diagU);
306:   if (nneg)  *nneg  = neg;
307:   if (nzero) *nzero = zero;
308:   if (npos)  *npos  = pos;
309:   return(0);
310: }

312: static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
313: {
314:   Mat_SuperLU_DIST  *lu = (Mat_SuperLU_DIST*)F->data;
315:   Mat               Aloc;
316:   const PetscScalar *av;
317:   const PetscInt    *ai=NULL,*aj=NULL;
318:   PetscInt          nz,dummy;
319:   int               sinfo;   /* SuperLU_Dist info flag is always an int even with long long indices */
320:   SuperLUStat_t     stat;
321:   double            *berr=0;
322:   PetscBool         ismpiaij,isseqaij,flg;
323:   PetscErrorCode    ierr;

326:   PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);
327:   PetscObjectBaseTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);
328:   if (ismpiaij) {
329:     MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&Aloc);
330:   } else if (isseqaij) {
331:     PetscObjectReference((PetscObject)A);
332:     Aloc = A;
333:   } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for type %s",((PetscObject)A)->type_name);

335:   MatGetRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
336:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GetRowIJ failed");
337:   MatSeqAIJGetArrayRead(Aloc,&av);
338:   nz   = ai[Aloc->rmap->n];

340:   /* Allocations for A_sup */
341:   if (lu->options.Fact == DOFACT) { /* first numeric factorization */
342: #if defined(PETSC_USE_COMPLEX)
343:     PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
344: #else
345:     PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
346: #endif
347:   } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
348:     if (lu->FactPattern == SamePattern_SameRowPerm) {
349:       lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
350:     } else if (lu->FactPattern == SamePattern) {
351:       PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct)); /* Deallocate L and U matrices. */
352:       lu->options.Fact = SamePattern;
353:     } else if (lu->FactPattern == DOFACT) {
354:       PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
355:       PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct));
356:       lu->options.Fact = DOFACT;

358: #if defined(PETSC_USE_COMPLEX)
359:       PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
360: #else
361:       PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
362: #endif
363:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm DOFACT");
364:   }

366:   /* Copy AIJ matrix to superlu_dist matrix */
367:   PetscArraycpy(lu->row,ai,Aloc->rmap->n+1);
368:   PetscArraycpy(lu->col,aj,nz);
369:   PetscArraycpy(lu->val,av,nz);
370:   MatRestoreRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
371:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"RestoreRowIJ failed");
372:   MatSeqAIJRestoreArrayRead(Aloc,&av);
373:   MatDestroy(&Aloc);

375:   /* Create and setup A_sup */
376:   if (lu->options.Fact == DOFACT) {
377: #if defined(PETSC_USE_COMPLEX)
378:     PetscStackCall("SuperLU_DIST:zCreate_CompRowLoc_Matrix_dist",zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, A->rmap->N, A->cmap->N, nz, A->rmap->n, A->rmap->rstart, lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE));
379: #else
380:     PetscStackCall("SuperLU_DIST:dCreate_CompRowLoc_Matrix_dist",dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, A->rmap->N, A->cmap->N, nz, A->rmap->n, A->rmap->rstart, lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE));
381: #endif
382:   }

384:   /* Factor the matrix. */
385:   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));   /* Initialize the statistics variables. */
386: #if defined(PETSC_USE_COMPLEX)
387:   PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, A->rmap->n, 0, &lu->grid, &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
388: #else
389:   PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, A->rmap->n, 0, &lu->grid, &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
390: #endif

392:   if (sinfo > 0) {
393:     if (A->erroriffailure) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
394:     else {
395:       if (sinfo <= lu->A_sup.ncol) {
396:         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
397:         PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);
398:       } else if (sinfo > lu->A_sup.ncol) {
399:         /*
400:          number of bytes allocated when memory allocation
401:          failure occurred, plus A->ncol.
402:          */
403:         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
404:         PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);
405:       }
406:     }
407:   } else if (sinfo < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo);

409:   if (lu->options.PrintStat) {
410:     PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid));  /* Print the statistics. */
411:   }
412:   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
413:   F->assembled     = PETSC_TRUE;
414:   F->preallocated  = PETSC_TRUE;
415:   lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
416:   return(0);
417: }

419: /* Note the Petsc r and c permutations are ignored */
420: static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
421: {
422:   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
423:   PetscInt         M   = A->rmap->N,N=A->cmap->N;

426:   /* Initialize ScalePermstruct and LUstruct. */
427:   PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
428:   PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
429:   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
430:   F->ops->solve           = MatSolve_SuperLU_DIST;
431:   F->ops->matsolve        = MatMatSolve_SuperLU_DIST;
432:   F->ops->getinertia      = NULL;

434:   if (A->symmetric || A->hermitian) F->ops->getinertia = MatGetInertia_SuperLU_DIST;
435:   lu->CleanUpSuperLU_Dist = PETSC_TRUE;
436:   return(0);
437: }

439: static PetscErrorCode MatCholeskyFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,const MatFactorInfo *info)
440: {

444:   MatLUFactorSymbolic_SuperLU_DIST(F,A,r,r,info);
445:   F->ops->choleskyfactornumeric = MatLUFactorNumeric_SuperLU_DIST;
446:   return(0);
447: }

449: static PetscErrorCode MatFactorGetSolverType_aij_superlu_dist(Mat A,MatSolverType *type)
450: {
452:   *type = MATSOLVERSUPERLU_DIST;
453:   return(0);
454: }

456: static PetscErrorCode MatView_Info_SuperLU_DIST(Mat A,PetscViewer viewer)
457: {
458:   Mat_SuperLU_DIST       *lu=(Mat_SuperLU_DIST*)A->data;
459:   superlu_dist_options_t options;
460:   PetscErrorCode         ierr;

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

466:   options = lu->options;
467:   PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");
468:   PetscViewerASCIIPrintf(viewer,"  Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);
469:   PetscViewerASCIIPrintf(viewer,"  Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);
470:   PetscViewerASCIIPrintf(viewer,"  Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);
471:   PetscViewerASCIIPrintf(viewer,"  Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);
472:   PetscViewerASCIIPrintf(viewer,"  Processors in row %d col partition %d \n",lu->nprow,lu->npcol);

474:   switch (options.RowPerm) {
475:   case NOROWPERM:
476:     PetscViewerASCIIPrintf(viewer,"  Row permutation NOROWPERM\n");
477:     break;
478:   case LargeDiag_MC64:
479:     PetscViewerASCIIPrintf(viewer,"  Row permutation LargeDiag_MC64\n");
480:     break;
481:   case LargeDiag_AWPM:
482:     PetscViewerASCIIPrintf(viewer,"  Row permutation LargeDiag_AWPM\n");
483:     break;
484:   case MY_PERMR:
485:     PetscViewerASCIIPrintf(viewer,"  Row permutation MY_PERMR\n");
486:     break;
487:   default:
488:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
489:   }

491:   switch (options.ColPerm) {
492:   case NATURAL:
493:     PetscViewerASCIIPrintf(viewer,"  Column permutation NATURAL\n");
494:     break;
495:   case MMD_AT_PLUS_A:
496:     PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_AT_PLUS_A\n");
497:     break;
498:   case MMD_ATA:
499:     PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_ATA\n");
500:     break;
501:   /*  Even though this is called METIS, the SuperLU_DIST code sets this by default if PARMETIS is defined, not METIS */
502:   case METIS_AT_PLUS_A:
503:     PetscViewerASCIIPrintf(viewer,"  Column permutation METIS_AT_PLUS_A\n");
504:     break;
505:   case PARMETIS:
506:     PetscViewerASCIIPrintf(viewer,"  Column permutation PARMETIS\n");
507:     break;
508:   default:
509:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
510:   }

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

514:   if (lu->FactPattern == SamePattern) {
515:     PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern\n");
516:   } else if (lu->FactPattern == SamePattern_SameRowPerm) {
517:     PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern_SameRowPerm\n");
518:   } else if (lu->FactPattern == DOFACT) {
519:     PetscViewerASCIIPrintf(viewer,"  Repeated factorization DOFACT\n");
520:   } else {
521:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown factorization pattern");
522:   }
523:   return(0);
524: }

526: static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
527: {
528:   PetscErrorCode    ierr;
529:   PetscBool         iascii;
530:   PetscViewerFormat format;

533:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
534:   if (iascii) {
535:     PetscViewerGetFormat(viewer,&format);
536:     if (format == PETSC_VIEWER_ASCII_INFO) {
537:       MatView_Info_SuperLU_DIST(A,viewer);
538:     }
539:   }
540:   return(0);
541: }

543: static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
544: {
545:   Mat                    B;
546:   Mat_SuperLU_DIST       *lu;
547:   PetscErrorCode         ierr;
548:   PetscInt               M=A->rmap->N,N=A->cmap->N,indx;
549:   PetscMPIInt            size;
550:   superlu_dist_options_t options;
551:   PetscBool              flg;
552:   const char             *colperm[]     = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
553:   const char             *rowperm[]     = {"NOROWPERM","LargeDiag_MC64","LargeDiag_AWPM","MY_PERMR"};
554:   const char             *factPattern[] = {"SamePattern","SamePattern_SameRowPerm","DOFACT"};
555:   PetscBool              set;

558:   /* Create the factorization matrix */
559:   MatCreate(PetscObjectComm((PetscObject)A),&B);
560:   MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);
561:   PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);
562:   MatSetUp(B);
563:   B->ops->getinfo = MatGetInfo_External;
564:   B->ops->view    = MatView_SuperLU_DIST;
565:   B->ops->destroy = MatDestroy_SuperLU_DIST;

567:   /* Set the default input options:
568:      options.Fact              = DOFACT;
569:      options.Equil             = YES;
570:      options.ParSymbFact       = NO;
571:      options.ColPerm           = METIS_AT_PLUS_A;
572:      options.RowPerm           = LargeDiag_MC64;
573:      options.ReplaceTinyPivot  = YES;
574:      options.IterRefine        = DOUBLE;
575:      options.Trans             = NOTRANS;
576:      options.SolveInitialized  = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
577:      options.RefineInitialized = NO;
578:      options.PrintStat         = YES;
579:      options.SymPattern        = NO;
580:   */
581:   set_default_options_dist(&options);

583:   if (ftype == MAT_FACTOR_LU) {
584:     B->factortype = MAT_FACTOR_LU;
585:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
586:   } else {
587:     B->factortype = MAT_FACTOR_CHOLESKY;
588:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SuperLU_DIST;
589:     options.SymPattern = YES;
590:   }

592:   /* set solvertype */
593:   PetscFree(B->solvertype);
594:   PetscStrallocpy(MATSOLVERSUPERLU_DIST,&B->solvertype);

596:   PetscNewLog(B,&lu);
597:   B->data = lu;
598:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);

600:   {
601:     PetscMPIInt       flg;
602:     MPI_Comm          comm;
603:     PetscSuperLU_DIST *context = NULL;

605:     PetscObjectGetComm((PetscObject)A,&comm);
606:     if (Petsc_Superlu_dist_keyval == MPI_KEYVAL_INVALID) {
607:       MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Superlu_dist_keyval_Delete_Fn,&Petsc_Superlu_dist_keyval,(void*)0);
608:       PetscRegisterFinalize(Petsc_Superlu_dist_keyval_free);
609:     }
610:     MPI_Comm_get_attr(comm,Petsc_Superlu_dist_keyval,&context,&flg);
611:     if (!flg || context->busy) {
612:       if (!flg) {
613:         PetscNew(&context);
614:         context->busy = PETSC_TRUE;
615:         MPI_Comm_dup(comm,&context->comm);
616:         MPI_Comm_set_attr(comm,Petsc_Superlu_dist_keyval,context);
617:       } else {
618:         MPI_Comm_dup(comm,&lu->comm_superlu);
619:       }

621:       /* Default num of process columns and rows */
622:       lu->nprow = (int_t) (0.5 + PetscSqrtReal((PetscReal)size));
623:       if (!lu->nprow) lu->nprow = 1;
624:       while (lu->nprow > 0) {
625:         lu->npcol = (int_t) (size/lu->nprow);
626:         if (size == lu->nprow * lu->npcol) break;
627:         lu->nprow--;
628:       }
629:       PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");
630:       PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);
631:       PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);
632:       PetscOptionsEnd();
633:       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);
634:       PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(context ? context->comm : lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid));
635:       if (context) context->grid = lu->grid;
636:       PetscInfo(NULL,"Duplicating a communicator for SuperLU_DIST and calling superlu_gridinit()\n");
637:       if (!flg) {
638:         PetscInfo(NULL,"Storing communicator and SuperLU_DIST grid in communicator attribute\n");
639:       } else {
640:         PetscInfo(NULL,"Communicator attribute already in use so not saving communicator and SuperLU_DIST grid in communicator attribute \n");
641:       }
642:     } else {
643:       PetscInfo(NULL,"Reusing communicator and superlu_gridinit() for SuperLU_DIST from communicator attribute.");
644:       context->busy = PETSC_TRUE;
645:       lu->grid      = context->grid;
646:     }
647:   }

649:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");
650:   PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
651:   if (set && !flg) options.Equil = NO;

653:   PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,4,rowperm[1],&indx,&flg);
654:   if (flg) {
655:     switch (indx) {
656:     case 0:
657:       options.RowPerm = NOROWPERM;
658:       break;
659:     case 1:
660:       options.RowPerm = LargeDiag_MC64;
661:       break;
662:     case 2:
663:       options.RowPerm = LargeDiag_AWPM;
664:       break;
665:     case 3:
666:       options.RowPerm = MY_PERMR;
667:       break;
668:     default:
669:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown row permutation");
670:     }
671:   }

673:   PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);
674:   if (flg) {
675:     switch (indx) {
676:     case 0:
677:       options.ColPerm = NATURAL;
678:       break;
679:     case 1:
680:       options.ColPerm = MMD_AT_PLUS_A;
681:       break;
682:     case 2:
683:       options.ColPerm = MMD_ATA;
684:       break;
685:     case 3:
686:       options.ColPerm = METIS_AT_PLUS_A;
687:       break;
688:     case 4:
689:       options.ColPerm = PARMETIS;   /* only works for np>1 */
690:       break;
691:     default:
692:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
693:     }
694:   }

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

700:   options.ParSymbFact = NO;
701:   PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);
702:   if (set && flg && size>1) {
703: #if defined(PETSC_HAVE_PARMETIS)
704:     options.ParSymbFact = YES;
705:     options.ColPerm     = PARMETIS;   /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */
706: #else
707:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS");
708: #endif
709:   }

711:   lu->FactPattern = SamePattern;
712:   PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,3,factPattern[0],&indx,&flg);
713:   if (flg) {
714:     switch (indx) {
715:     case 0:
716:       lu->FactPattern = SamePattern;
717:       break;
718:     case 1:
719:       lu->FactPattern = SamePattern_SameRowPerm;
720:       break;
721:     case 2:
722:       lu->FactPattern = DOFACT;
723:       break;
724:     }
725:   }

727:   options.IterRefine = NOREFINE;
728:   PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);
729:   if (set) {
730:     if (flg) options.IterRefine = SLU_DOUBLE;
731:     else options.IterRefine = NOREFINE;
732:   }

734:   if (PetscLogPrintInfo) options.PrintStat = YES;
735:   else options.PrintStat = NO;
736:   PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);
737:   PetscOptionsEnd();

739:   lu->options              = options;
740:   lu->options.Fact         = DOFACT;
741:   lu->matsolve_iscalled    = PETSC_FALSE;
742:   lu->matmatsolve_iscalled = PETSC_FALSE;

744:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_superlu_dist);
745:   PetscObjectComposeFunction((PetscObject)B,"MatSuperluDistGetDiagU_C",MatSuperluDistGetDiagU_SuperLU_DIST);

747:   *F = B;
748:   return(0);
749: }

751: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU_DIST(void)
752: {
755:   MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
756:   MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
757:   MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);
758:   MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);
759:   return(0);
760: }

762: /*MC
763:   MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization

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

767:   Use -pc_type lu -pc_factor_mat_solver_type superlu_dist to use this direct solver

769:    Works with AIJ matrices

771:   Options Database Keys:
772: + -mat_superlu_dist_r <n> - number of rows in processor partition
773: . -mat_superlu_dist_c <n> - number of columns in processor partition
774: . -mat_superlu_dist_equil - equilibrate the matrix
775: . -mat_superlu_dist_rowperm <NOROWPERM,LargeDiag_MC64,LargeDiag_AWPM,MY_PERMR> - row permutation
776: . -mat_superlu_dist_colperm <NATURAL,MMD_AT_PLUS_A,MMD_ATA,METIS_AT_PLUS_A,PARMETIS> - column permutation
777: . -mat_superlu_dist_replacetinypivot - replace tiny pivots
778: . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm DOFACT
779: . -mat_superlu_dist_iterrefine - use iterative refinement
780: - -mat_superlu_dist_statprint - print factorization information

782:    Level: beginner

784: .seealso: PCLU

786: .seealso: PCFactorSetMatSolverType(), MatSolverType

788: M*/