Actual source code: superlu_dist.c
petsc-3.14.2 2020-12-03
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: static PetscErrorCode Petsc_Superlu_dist_keyval_free(void)
104: {
108: PetscInfo(NULL,"Freeing Petsc_Superlu_dist_keyval\n");
109: MPI_Comm_free_keyval(&Petsc_Superlu_dist_keyval);
110: return(0);
111: }
113: static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
114: {
115: PetscErrorCode ierr;
116: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
119: if (lu->CleanUpSuperLU_Dist) {
120: /* Deallocate SuperLU_DIST storage */
121: PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
122: if (lu->options.SolveInitialized) {
123: #if defined(PETSC_USE_COMPLEX)
124: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
125: #else
126: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
127: #endif
128: }
129: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
130: PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
131: PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));
133: /* Release the SuperLU_DIST process grid. Only if the matrix has its own copy, this is it is not in the communicator context */
134: if (lu->comm_superlu) {
135: PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
136: MPI_Comm_free(&(lu->comm_superlu));
137: } else {
138: PetscSuperLU_DIST *context;
139: MPI_Comm comm;
140: PetscMPIInt flg;
142: PetscObjectGetComm((PetscObject)A,&comm);
143: MPI_Comm_get_attr(comm,Petsc_Superlu_dist_keyval,&context,&flg);
144: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Communicator does not have expected Petsc_Superlu_dist_keyval attribute");
145: context->busy = PETSC_FALSE;
146: }
147: }
148: PetscFree(A->data);
149: /* clear composed functions */
150: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
151: PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);
153: return(0);
154: }
156: static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
157: {
158: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
159: PetscErrorCode ierr;
160: PetscInt m=A->rmap->n;
161: SuperLUStat_t stat;
162: double berr[1];
163: PetscScalar *bptr=NULL;
164: int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
165: static PetscBool cite = PETSC_FALSE;
168: if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
169: 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);
171: if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
172: /* see comments in MatMatSolve() */
173: #if defined(PETSC_USE_COMPLEX)
174: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
175: #else
176: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
177: #endif
178: lu->options.SolveInitialized = NO;
179: }
180: VecCopy(b_mpi,x);
181: VecGetArray(x,&bptr);
183: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
184: #if defined(PETSC_USE_COMPLEX)
185: 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));
186: #else
187: PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,1,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
188: #endif
189: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
191: if (lu->options.PrintStat) PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
192: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
194: VecRestoreArray(x,&bptr);
195: lu->matsolve_iscalled = PETSC_TRUE;
196: lu->matmatsolve_iscalled = PETSC_FALSE;
197: return(0);
198: }
200: static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
201: {
202: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
203: PetscErrorCode ierr;
204: PetscInt m=A->rmap->n,nrhs;
205: SuperLUStat_t stat;
206: double berr[1];
207: PetscScalar *bptr;
208: int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
209: PetscBool flg;
212: if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
213: PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
214: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
215: if (X != B_mpi) {
216: PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
217: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
218: }
220: if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
221: /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
222: thus destroy it and create a new SOLVEstruct.
223: Otherwise it may result in memory corruption or incorrect solution
224: See src/mat/tests/ex125.c */
225: #if defined(PETSC_USE_COMPLEX)
226: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
227: #else
228: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
229: #endif
230: lu->options.SolveInitialized = NO;
231: }
232: if (X != B_mpi) {
233: MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);
234: }
236: MatGetSize(B_mpi,NULL,&nrhs);
238: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
239: MatDenseGetArray(X,&bptr);
241: #if defined(PETSC_USE_COMPLEX)
242: 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));
243: #else
244: PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
245: #endif
247: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
248: MatDenseRestoreArray(X,&bptr);
250: if (lu->options.PrintStat) PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
251: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
252: lu->matsolve_iscalled = PETSC_FALSE;
253: lu->matmatsolve_iscalled = PETSC_TRUE;
254: return(0);
255: }
257: /*
258: input:
259: F: numeric Cholesky factor
260: output:
261: nneg: total number of negative pivots
262: nzero: total number of zero pivots
263: npos: (global dimension of F) - nneg - nzero
264: */
265: static PetscErrorCode MatGetInertia_SuperLU_DIST(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
266: {
267: PetscErrorCode ierr;
268: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
269: PetscScalar *diagU=NULL;
270: PetscInt M,i,neg=0,zero=0,pos=0;
271: PetscReal r;
274: if (!F->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix factor F is not assembled");
275: if (lu->options.RowPerm != NOROWPERM) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must set NOROWPERM");
276: MatGetSize(F,&M,NULL);
277: PetscMalloc1(M,&diagU);
278: MatSuperluDistGetDiagU(F,diagU);
279: for (i=0; i<M; i++) {
280: #if defined(PETSC_USE_COMPLEX)
281: r = PetscImaginaryPart(diagU[i])/10.0;
282: 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);
283: r = PetscRealPart(diagU[i]);
284: #else
285: r = diagU[i];
286: #endif
287: if (r > 0) {
288: pos++;
289: } else if (r < 0) {
290: neg++;
291: } else zero++;
292: }
294: PetscFree(diagU);
295: if (nneg) *nneg = neg;
296: if (nzero) *nzero = zero;
297: if (npos) *npos = pos;
298: return(0);
299: }
301: static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
302: {
303: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
304: Mat Aloc;
305: const PetscScalar *av;
306: const PetscInt *ai=NULL,*aj=NULL;
307: PetscInt nz,dummy;
308: int sinfo; /* SuperLU_Dist info flag is always an int even with long long indices */
309: SuperLUStat_t stat;
310: double *berr=0;
311: PetscBool ismpiaij,isseqaij,flg;
312: PetscErrorCode ierr;
315: PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);
316: PetscObjectBaseTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);
317: if (ismpiaij) {
318: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&Aloc);
319: } else if (isseqaij) {
320: PetscObjectReference((PetscObject)A);
321: Aloc = A;
322: } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for type %s",((PetscObject)A)->type_name);
324: MatGetRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
325: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GetRowIJ failed");
326: MatSeqAIJGetArrayRead(Aloc,&av);
327: nz = ai[Aloc->rmap->n];
329: /* Allocations for A_sup */
330: if (lu->options.Fact == DOFACT) { /* first numeric factorization */
331: #if defined(PETSC_USE_COMPLEX)
332: PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
333: #else
334: PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
335: #endif
336: } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
337: if (lu->FactPattern == SamePattern_SameRowPerm) {
338: lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
339: } else if (lu->FactPattern == SamePattern) {
340: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct)); /* Deallocate L and U matrices. */
341: lu->options.Fact = SamePattern;
342: } else if (lu->FactPattern == DOFACT) {
343: PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
344: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct));
345: lu->options.Fact = DOFACT;
347: #if defined(PETSC_USE_COMPLEX)
348: PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
349: #else
350: PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
351: #endif
352: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm DOFACT");
353: }
355: /* Copy AIJ matrix to superlu_dist matrix */
356: PetscArraycpy(lu->row,ai,Aloc->rmap->n+1);
357: PetscArraycpy(lu->col,aj,nz);
358: PetscArraycpy(lu->val,av,nz);
359: MatRestoreRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
360: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"RestoreRowIJ failed");
361: MatSeqAIJRestoreArrayRead(Aloc,&av);
362: MatDestroy(&Aloc);
364: /* Create and setup A_sup */
365: if (lu->options.Fact == DOFACT) {
366: #if defined(PETSC_USE_COMPLEX)
367: 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));
368: #else
369: 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));
370: #endif
371: }
373: /* Factor the matrix. */
374: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
375: #if defined(PETSC_USE_COMPLEX)
376: 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));
377: #else
378: 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));
379: #endif
381: if (sinfo > 0) {
382: if (A->erroriffailure) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
383: else {
384: if (sinfo <= lu->A_sup.ncol) {
385: F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
386: PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);
387: } else if (sinfo > lu->A_sup.ncol) {
388: /*
389: number of bytes allocated when memory allocation
390: failure occurred, plus A->ncol.
391: */
392: F->factorerrortype = MAT_FACTOR_OUTMEMORY;
393: PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);
394: }
395: }
396: } else if (sinfo < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo);
398: if (lu->options.PrintStat) {
399: PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
400: }
401: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
402: F->assembled = PETSC_TRUE;
403: F->preallocated = PETSC_TRUE;
404: lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
405: return(0);
406: }
408: /* Note the Petsc r and c permutations are ignored */
409: static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
410: {
411: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
412: PetscInt M = A->rmap->N,N=A->cmap->N;
415: /* Initialize ScalePermstruct and LUstruct. */
416: PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
417: PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
418: F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
419: F->ops->solve = MatSolve_SuperLU_DIST;
420: F->ops->matsolve = MatMatSolve_SuperLU_DIST;
421: F->ops->getinertia = NULL;
423: if (A->symmetric || A->hermitian) F->ops->getinertia = MatGetInertia_SuperLU_DIST;
424: lu->CleanUpSuperLU_Dist = PETSC_TRUE;
425: return(0);
426: }
428: static PetscErrorCode MatCholeskyFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,const MatFactorInfo *info)
429: {
433: MatLUFactorSymbolic_SuperLU_DIST(F,A,r,r,info);
434: F->ops->choleskyfactornumeric = MatLUFactorNumeric_SuperLU_DIST;
435: return(0);
436: }
438: static PetscErrorCode MatFactorGetSolverType_aij_superlu_dist(Mat A,MatSolverType *type)
439: {
441: *type = MATSOLVERSUPERLU_DIST;
442: return(0);
443: }
445: static PetscErrorCode MatView_Info_SuperLU_DIST(Mat A,PetscViewer viewer)
446: {
447: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)A->data;
448: superlu_dist_options_t options;
449: PetscErrorCode ierr;
452: /* check if matrix is superlu_dist type */
453: if (A->ops->solve != MatSolve_SuperLU_DIST) return(0);
455: options = lu->options;
456: PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");
457: PetscViewerASCIIPrintf(viewer," Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);
458: PetscViewerASCIIPrintf(viewer," Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);
459: PetscViewerASCIIPrintf(viewer," Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);
460: PetscViewerASCIIPrintf(viewer," Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);
461: PetscViewerASCIIPrintf(viewer," Processors in row %d col partition %d \n",lu->nprow,lu->npcol);
463: switch (options.RowPerm) {
464: case NOROWPERM:
465: PetscViewerASCIIPrintf(viewer," Row permutation NOROWPERM\n");
466: break;
467: case LargeDiag_MC64:
468: PetscViewerASCIIPrintf(viewer," Row permutation LargeDiag_MC64\n");
469: break;
470: case LargeDiag_AWPM:
471: PetscViewerASCIIPrintf(viewer," Row permutation LargeDiag_AWPM\n");
472: break;
473: case MY_PERMR:
474: PetscViewerASCIIPrintf(viewer," Row permutation MY_PERMR\n");
475: break;
476: default:
477: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
478: }
480: switch (options.ColPerm) {
481: case NATURAL:
482: PetscViewerASCIIPrintf(viewer," Column permutation NATURAL\n");
483: break;
484: case MMD_AT_PLUS_A:
485: PetscViewerASCIIPrintf(viewer," Column permutation MMD_AT_PLUS_A\n");
486: break;
487: case MMD_ATA:
488: PetscViewerASCIIPrintf(viewer," Column permutation MMD_ATA\n");
489: break;
490: /* Even though this is called METIS, the SuperLU_DIST code sets this by default if PARMETIS is defined, not METIS */
491: case METIS_AT_PLUS_A:
492: PetscViewerASCIIPrintf(viewer," Column permutation METIS_AT_PLUS_A\n");
493: break;
494: case PARMETIS:
495: PetscViewerASCIIPrintf(viewer," Column permutation PARMETIS\n");
496: break;
497: default:
498: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
499: }
501: PetscViewerASCIIPrintf(viewer," Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);
503: if (lu->FactPattern == SamePattern) {
504: PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern\n");
505: } else if (lu->FactPattern == SamePattern_SameRowPerm) {
506: PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern_SameRowPerm\n");
507: } else if (lu->FactPattern == DOFACT) {
508: PetscViewerASCIIPrintf(viewer," Repeated factorization DOFACT\n");
509: } else {
510: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown factorization pattern");
511: }
512: return(0);
513: }
515: static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
516: {
517: PetscErrorCode ierr;
518: PetscBool iascii;
519: PetscViewerFormat format;
522: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
523: if (iascii) {
524: PetscViewerGetFormat(viewer,&format);
525: if (format == PETSC_VIEWER_ASCII_INFO) {
526: MatView_Info_SuperLU_DIST(A,viewer);
527: }
528: }
529: return(0);
530: }
532: static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
533: {
534: Mat B;
535: Mat_SuperLU_DIST *lu;
536: PetscErrorCode ierr;
537: PetscInt M=A->rmap->N,N=A->cmap->N,indx;
538: PetscMPIInt size;
539: superlu_dist_options_t options;
540: PetscBool flg;
541: const char *colperm[] = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
542: const char *rowperm[] = {"NOROWPERM","LargeDiag_MC64","LargeDiag_AWPM","MY_PERMR"};
543: const char *factPattern[] = {"SamePattern","SamePattern_SameRowPerm","DOFACT"};
544: PetscBool set;
547: /* Create the factorization matrix */
548: MatCreate(PetscObjectComm((PetscObject)A),&B);
549: MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);
550: PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);
551: MatSetUp(B);
552: B->ops->getinfo = MatGetInfo_External;
553: B->ops->view = MatView_SuperLU_DIST;
554: B->ops->destroy = MatDestroy_SuperLU_DIST;
556: /* Set the default input options:
557: options.Fact = DOFACT;
558: options.Equil = YES;
559: options.ParSymbFact = NO;
560: options.ColPerm = METIS_AT_PLUS_A;
561: options.RowPerm = LargeDiag_MC64;
562: options.ReplaceTinyPivot = YES;
563: options.IterRefine = DOUBLE;
564: options.Trans = NOTRANS;
565: options.SolveInitialized = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
566: options.RefineInitialized = NO;
567: options.PrintStat = YES;
568: options.SymPattern = NO;
569: */
570: set_default_options_dist(&options);
572: if (ftype == MAT_FACTOR_LU) {
573: B->factortype = MAT_FACTOR_LU;
574: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
575: } else {
576: B->factortype = MAT_FACTOR_CHOLESKY;
577: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SuperLU_DIST;
578: options.SymPattern = YES;
579: }
581: /* set solvertype */
582: PetscFree(B->solvertype);
583: PetscStrallocpy(MATSOLVERSUPERLU_DIST,&B->solvertype);
585: PetscNewLog(B,&lu);
586: B->data = lu;
587: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
589: {
590: PetscMPIInt flg;
591: MPI_Comm comm;
592: PetscSuperLU_DIST *context = NULL;
594: PetscObjectGetComm((PetscObject)A,&comm);
595: if (Petsc_Superlu_dist_keyval == MPI_KEYVAL_INVALID) {
596: MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Superlu_dist_keyval_Delete_Fn,&Petsc_Superlu_dist_keyval,(void*)0);
597: PetscRegisterFinalize(Petsc_Superlu_dist_keyval_free);
598: }
599: MPI_Comm_get_attr(comm,Petsc_Superlu_dist_keyval,&context,&flg);
600: if (!flg || context->busy) {
601: if (!flg) {
602: PetscNew(&context);
603: context->busy = PETSC_TRUE;
604: MPI_Comm_dup(comm,&context->comm);
605: MPI_Comm_set_attr(comm,Petsc_Superlu_dist_keyval,context);
606: } else {
607: MPI_Comm_dup(comm,&lu->comm_superlu);
608: }
610: /* Default num of process columns and rows */
611: lu->nprow = (int_t) (0.5 + PetscSqrtReal((PetscReal)size));
612: if (!lu->nprow) lu->nprow = 1;
613: while (lu->nprow > 0) {
614: lu->npcol = (int_t) (size/lu->nprow);
615: if (size == lu->nprow * lu->npcol) break;
616: lu->nprow--;
617: }
618: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");
619: PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);
620: PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);
621: PetscOptionsEnd();
622: 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);
623: PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(context ? context->comm : lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid));
624: if (context) context->grid = lu->grid;
625: PetscInfo(NULL,"Duplicating a communicator for SuperLU_DIST and calling superlu_gridinit()\n");
626: if (!flg) {
627: PetscInfo(NULL,"Storing communicator and SuperLU_DIST grid in communicator attribute\n");
628: } else {
629: PetscInfo(NULL,"Communicator attribute already in use so not saving communicator and SuperLU_DIST grid in communicator attribute \n");
630: }
631: } else {
632: PetscInfo(NULL,"Reusing communicator and superlu_gridinit() for SuperLU_DIST from communicator attribute.");
633: context->busy = PETSC_TRUE;
634: lu->grid = context->grid;
635: }
636: }
638: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");
639: PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
640: if (set && !flg) options.Equil = NO;
642: PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,4,rowperm[1],&indx,&flg);
643: if (flg) {
644: switch (indx) {
645: case 0:
646: options.RowPerm = NOROWPERM;
647: break;
648: case 1:
649: options.RowPerm = LargeDiag_MC64;
650: break;
651: case 2:
652: options.RowPerm = LargeDiag_AWPM;
653: break;
654: case 3:
655: options.RowPerm = MY_PERMR;
656: break;
657: default:
658: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown row permutation");
659: }
660: }
662: PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);
663: if (flg) {
664: switch (indx) {
665: case 0:
666: options.ColPerm = NATURAL;
667: break;
668: case 1:
669: options.ColPerm = MMD_AT_PLUS_A;
670: break;
671: case 2:
672: options.ColPerm = MMD_ATA;
673: break;
674: case 3:
675: options.ColPerm = METIS_AT_PLUS_A;
676: break;
677: case 4:
678: options.ColPerm = PARMETIS; /* only works for np>1 */
679: break;
680: default:
681: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
682: }
683: }
685: options.ReplaceTinyPivot = NO;
686: PetscOptionsBool("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",options.ReplaceTinyPivot ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
687: if (set && flg) options.ReplaceTinyPivot = YES;
689: options.ParSymbFact = NO;
690: PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);
691: if (set && flg && size>1) {
692: #if defined(PETSC_HAVE_PARMETIS)
693: options.ParSymbFact = YES;
694: options.ColPerm = PARMETIS; /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */
695: #else
696: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS");
697: #endif
698: }
700: lu->FactPattern = SamePattern;
701: PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,3,factPattern[0],&indx,&flg);
702: if (flg) {
703: switch (indx) {
704: case 0:
705: lu->FactPattern = SamePattern;
706: break;
707: case 1:
708: lu->FactPattern = SamePattern_SameRowPerm;
709: break;
710: case 2:
711: lu->FactPattern = DOFACT;
712: break;
713: }
714: }
716: options.IterRefine = NOREFINE;
717: PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);
718: if (set) {
719: if (flg) options.IterRefine = SLU_DOUBLE;
720: else options.IterRefine = NOREFINE;
721: }
723: if (PetscLogPrintInfo) options.PrintStat = YES;
724: else options.PrintStat = NO;
725: PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);
726: PetscOptionsEnd();
728: lu->options = options;
729: lu->options.Fact = DOFACT;
730: lu->matsolve_iscalled = PETSC_FALSE;
731: lu->matmatsolve_iscalled = PETSC_FALSE;
733: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_superlu_dist);
734: PetscObjectComposeFunction((PetscObject)B,"MatSuperluDistGetDiagU_C",MatSuperluDistGetDiagU_SuperLU_DIST);
736: *F = B;
737: return(0);
738: }
740: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU_DIST(void)
741: {
744: MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
745: MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);
746: MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);
747: MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);
748: return(0);
749: }
751: /*MC
752: MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization
754: Use ./configure --download-superlu_dist --download-parmetis --download-metis --download-ptscotch to have PETSc installed with SuperLU_DIST
756: Use -pc_type lu -pc_factor_mat_solver_type superlu_dist to use this direct solver
758: Works with AIJ matrices
760: Options Database Keys:
761: + -mat_superlu_dist_r <n> - number of rows in processor partition
762: . -mat_superlu_dist_c <n> - number of columns in processor partition
763: . -mat_superlu_dist_equil - equilibrate the matrix
764: . -mat_superlu_dist_rowperm <NOROWPERM,LargeDiag_MC64,LargeDiag_AWPM,MY_PERMR> - row permutation
765: . -mat_superlu_dist_colperm <NATURAL,MMD_AT_PLUS_A,MMD_ATA,METIS_AT_PLUS_A,PARMETIS> - column permutation
766: . -mat_superlu_dist_replacetinypivot - replace tiny pivots
767: . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm DOFACT
768: . -mat_superlu_dist_iterrefine - use iterative refinement
769: - -mat_superlu_dist_statprint - print factorization information
771: Level: beginner
773: .seealso: PCLU
775: .seealso: PCFactorSetMatSolverType(), MatSolverType
777: M*/