Actual source code: superlu_dist.c
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;
57: PetscErrorCode MatSuperluDistGetDiagU_SuperLU_DIST(Mat F,PetscScalar *diagU)
58: {
59: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
62: #if defined(PETSC_USE_COMPLEX)
63: PetscStackCall("SuperLU_DIST:pzGetDiagU",pzGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,(doublecomplex*)diagU));
64: #else
65: PetscStackCall("SuperLU_DIST:pdGetDiagU",pdGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,diagU));
66: #endif
67: return(0);
68: }
70: PetscErrorCode MatSuperluDistGetDiagU(Mat F,PetscScalar *diagU)
71: {
76: PetscTryMethod(F,"MatSuperluDistGetDiagU_C",(Mat,PetscScalar*),(F,diagU));
77: return(0);
78: }
80: /* This allows reusing the Superlu_DIST communicator and grid when only a single SuperLU_DIST matrix is used at a time */
81: typedef struct {
82: MPI_Comm comm;
83: PetscBool busy;
84: gridinfo_t grid;
85: } PetscSuperLU_DIST;
86: static PetscMPIInt Petsc_Superlu_dist_keyval = MPI_KEYVAL_INVALID;
88: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Superlu_dist_keyval_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
89: {
90: PetscErrorCode ierr;
91: PetscSuperLU_DIST *context = (PetscSuperLU_DIST *) attr_val;
94: if (keyval != Petsc_Superlu_dist_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
95: PetscInfo(NULL,"Removing Petsc_Superlu_dist_keyval attribute from communicator that is being freed\n");
96: PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&context->grid));
97: MPI_Comm_free(&context->comm);
98: PetscFree(context);
99: PetscFunctionReturn(MPI_SUCCESS);
100: }
102: /*
103: Performs MPI_Comm_free_keyval() on Petsc_Superlu_dist_keyval but keeps the global variable for
104: users who do not destroy all PETSc objects before PetscFinalize().
106: The value Petsc_Superlu_dist_keyval is retained so that Petsc_Superlu_dist_keyval_Delete_Fn()
107: can still check that the keyval associated with the MPI communicator is correct when the MPI
108: communicator is destroyed.
110: This is called in PetscFinalize()
111: */
112: static PetscErrorCode Petsc_Superlu_dist_keyval_free(void)
113: {
115: PetscMPIInt Petsc_Superlu_dist_keyval_temp = Petsc_Superlu_dist_keyval;
118: PetscInfo(NULL,"Freeing Petsc_Superlu_dist_keyval\n");
119: MPI_Comm_free_keyval(&Petsc_Superlu_dist_keyval_temp);
120: return(0);
121: }
123: static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
124: {
125: PetscErrorCode ierr;
126: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
129: if (lu->CleanUpSuperLU_Dist) {
130: /* Deallocate SuperLU_DIST storage */
131: PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
132: if (lu->options.SolveInitialized) {
133: #if defined(PETSC_USE_COMPLEX)
134: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
135: #else
136: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
137: #endif
138: }
139: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
140: PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
141: PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));
143: /* Release the SuperLU_DIST process grid. Only if the matrix has its own copy, this is it is not in the communicator context */
144: if (lu->comm_superlu) {
145: PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
146: MPI_Comm_free(&(lu->comm_superlu));
147: } else {
148: PetscSuperLU_DIST *context;
149: MPI_Comm comm;
150: PetscMPIInt flg;
152: PetscObjectGetComm((PetscObject)A,&comm);
153: MPI_Comm_get_attr(comm,Petsc_Superlu_dist_keyval,&context,&flg);
154: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Communicator does not have expected Petsc_Superlu_dist_keyval attribute");
155: context->busy = PETSC_FALSE;
156: }
157: }
158: PetscFree(A->data);
159: /* clear composed functions */
160: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
161: PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);
163: return(0);
164: }
166: static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
167: {
168: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
169: PetscErrorCode ierr;
170: PetscInt m=A->rmap->n;
171: SuperLUStat_t stat;
172: double berr[1];
173: PetscScalar *bptr=NULL;
174: int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
175: static PetscBool cite = PETSC_FALSE;
178: if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
179: 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);
181: if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
182: /* see comments in MatMatSolve() */
183: #if defined(PETSC_USE_COMPLEX)
184: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
185: #else
186: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
187: #endif
188: lu->options.SolveInitialized = NO;
189: }
190: VecCopy(b_mpi,x);
191: VecGetArray(x,&bptr);
193: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
194: #if defined(PETSC_USE_COMPLEX)
195: 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));
196: #else
197: PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,1,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
198: #endif
199: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
201: if (lu->options.PrintStat) PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
202: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
204: VecRestoreArray(x,&bptr);
205: lu->matsolve_iscalled = PETSC_TRUE;
206: lu->matmatsolve_iscalled = PETSC_FALSE;
207: return(0);
208: }
210: static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
211: {
212: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
213: PetscErrorCode ierr;
214: PetscInt m=A->rmap->n,nrhs;
215: SuperLUStat_t stat;
216: double berr[1];
217: PetscScalar *bptr;
218: int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
219: PetscBool flg;
222: if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact must equal FACTORED");
223: PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
224: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
225: if (X != B_mpi) {
226: PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
227: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
228: }
230: if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
231: /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
232: thus destroy it and create a new SOLVEstruct.
233: Otherwise it may result in memory corruption or incorrect solution
234: See src/mat/tests/ex125.c */
235: #if defined(PETSC_USE_COMPLEX)
236: PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
237: #else
238: PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
239: #endif
240: lu->options.SolveInitialized = NO;
241: }
242: if (X != B_mpi) {
243: MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);
244: }
246: MatGetSize(B_mpi,NULL,&nrhs);
248: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
249: MatDenseGetArray(X,&bptr);
251: #if defined(PETSC_USE_COMPLEX)
252: 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));
253: #else
254: PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
255: #endif
257: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
258: MatDenseRestoreArray(X,&bptr);
260: if (lu->options.PrintStat) PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
261: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
262: lu->matsolve_iscalled = PETSC_FALSE;
263: lu->matmatsolve_iscalled = PETSC_TRUE;
264: return(0);
265: }
267: /*
268: input:
269: F: numeric Cholesky factor
270: output:
271: nneg: total number of negative pivots
272: nzero: total number of zero pivots
273: npos: (global dimension of F) - nneg - nzero
274: */
275: static PetscErrorCode MatGetInertia_SuperLU_DIST(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
276: {
277: PetscErrorCode ierr;
278: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
279: PetscScalar *diagU=NULL;
280: PetscInt M,i,neg=0,zero=0,pos=0;
281: PetscReal r;
284: if (!F->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix factor F is not assembled");
285: if (lu->options.RowPerm != NOROWPERM) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must set NOROWPERM");
286: MatGetSize(F,&M,NULL);
287: PetscMalloc1(M,&diagU);
288: MatSuperluDistGetDiagU(F,diagU);
289: for (i=0; i<M; i++) {
290: #if defined(PETSC_USE_COMPLEX)
291: r = PetscImaginaryPart(diagU[i])/10.0;
292: 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);
293: r = PetscRealPart(diagU[i]);
294: #else
295: r = diagU[i];
296: #endif
297: if (r > 0) {
298: pos++;
299: } else if (r < 0) {
300: neg++;
301: } else zero++;
302: }
304: PetscFree(diagU);
305: if (nneg) *nneg = neg;
306: if (nzero) *nzero = zero;
307: if (npos) *npos = pos;
308: return(0);
309: }
311: static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
312: {
313: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
314: Mat Aloc;
315: const PetscScalar *av;
316: const PetscInt *ai=NULL,*aj=NULL;
317: PetscInt nz,dummy;
318: int sinfo; /* SuperLU_Dist info flag is always an int even with long long indices */
319: SuperLUStat_t stat;
320: double *berr=0;
321: PetscBool ismpiaij,isseqaij,flg;
322: PetscErrorCode ierr;
325: PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);
326: PetscObjectBaseTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);
327: if (ismpiaij) {
328: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&Aloc);
329: } else if (isseqaij) {
330: PetscObjectReference((PetscObject)A);
331: Aloc = A;
332: } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for type %s",((PetscObject)A)->type_name);
334: MatGetRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
335: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GetRowIJ failed");
336: MatSeqAIJGetArrayRead(Aloc,&av);
337: nz = ai[Aloc->rmap->n];
339: /* Allocations for A_sup */
340: if (lu->options.Fact == DOFACT) { /* first numeric factorization */
341: #if defined(PETSC_USE_COMPLEX)
342: PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
343: #else
344: PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
345: #endif
346: } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
347: if (lu->FactPattern == SamePattern_SameRowPerm) {
348: lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
349: } else if (lu->FactPattern == SamePattern) {
350: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct)); /* Deallocate L and U matrices. */
351: lu->options.Fact = SamePattern;
352: } else if (lu->FactPattern == DOFACT) {
353: PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
354: PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->rmap->N, &lu->grid, &lu->LUstruct));
355: lu->options.Fact = DOFACT;
357: #if defined(PETSC_USE_COMPLEX)
358: PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
359: #else
360: PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(Aloc->rmap->n, nz, &lu->val, &lu->col, &lu->row));
361: #endif
362: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm DOFACT");
363: }
365: /* Copy AIJ matrix to superlu_dist matrix */
366: PetscArraycpy(lu->row,ai,Aloc->rmap->n+1);
367: PetscArraycpy(lu->col,aj,nz);
368: PetscArraycpy(lu->val,av,nz);
369: MatRestoreRowIJ(Aloc,0,PETSC_FALSE,PETSC_FALSE,&dummy,&ai,&aj,&flg);
370: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"RestoreRowIJ failed");
371: MatSeqAIJRestoreArrayRead(Aloc,&av);
372: MatDestroy(&Aloc);
374: /* Create and setup A_sup */
375: if (lu->options.Fact == DOFACT) {
376: #if defined(PETSC_USE_COMPLEX)
377: 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));
378: #else
379: 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));
380: #endif
381: }
383: /* Factor the matrix. */
384: PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */
385: #if defined(PETSC_USE_COMPLEX)
386: 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));
387: #else
388: 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));
389: #endif
391: if (sinfo > 0) {
392: if (A->erroriffailure) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
393: else {
394: if (sinfo <= lu->A_sup.ncol) {
395: F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
396: PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);
397: } else if (sinfo > lu->A_sup.ncol) {
398: /*
399: number of bytes allocated when memory allocation
400: failure occurred, plus A->ncol.
401: */
402: F->factorerrortype = MAT_FACTOR_OUTMEMORY;
403: PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);
404: }
405: }
406: } else if (sinfo < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo);
408: if (lu->options.PrintStat) {
409: PetscStackCall("SuperLU_DIST:PStatPrint",PStatPrint(&lu->options, &stat, &lu->grid)); /* Print the statistics. */
410: }
411: PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
412: F->assembled = PETSC_TRUE;
413: F->preallocated = PETSC_TRUE;
414: lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
415: return(0);
416: }
418: /* Note the Petsc r and c permutations are ignored */
419: static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
420: {
421: Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
422: PetscInt M = A->rmap->N,N=A->cmap->N;
425: /* Initialize ScalePermstruct and LUstruct. */
426: PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
427: PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
428: F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
429: F->ops->solve = MatSolve_SuperLU_DIST;
430: F->ops->matsolve = MatMatSolve_SuperLU_DIST;
431: F->ops->getinertia = NULL;
433: if (A->symmetric || A->hermitian) F->ops->getinertia = MatGetInertia_SuperLU_DIST;
434: lu->CleanUpSuperLU_Dist = PETSC_TRUE;
435: return(0);
436: }
438: static PetscErrorCode MatCholeskyFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,const MatFactorInfo *info)
439: {
443: MatLUFactorSymbolic_SuperLU_DIST(F,A,r,r,info);
444: F->ops->choleskyfactornumeric = MatLUFactorNumeric_SuperLU_DIST;
445: return(0);
446: }
448: static PetscErrorCode MatFactorGetSolverType_aij_superlu_dist(Mat A,MatSolverType *type)
449: {
451: *type = MATSOLVERSUPERLU_DIST;
452: return(0);
453: }
455: static PetscErrorCode MatView_Info_SuperLU_DIST(Mat A,PetscViewer viewer)
456: {
457: Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)A->data;
458: superlu_dist_options_t options;
459: PetscErrorCode ierr;
462: /* check if matrix is superlu_dist type */
463: if (A->ops->solve != MatSolve_SuperLU_DIST) return(0);
465: options = lu->options;
466: PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");
467: PetscViewerASCIIPrintf(viewer," Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);
468: PetscViewerASCIIPrintf(viewer," Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);
469: PetscViewerASCIIPrintf(viewer," Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);
470: PetscViewerASCIIPrintf(viewer," Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);
471: PetscViewerASCIIPrintf(viewer," Processors in row %d col partition %d \n",lu->nprow,lu->npcol);
473: switch (options.RowPerm) {
474: case NOROWPERM:
475: PetscViewerASCIIPrintf(viewer," Row permutation NOROWPERM\n");
476: break;
477: case LargeDiag_MC64:
478: PetscViewerASCIIPrintf(viewer," Row permutation LargeDiag_MC64\n");
479: break;
480: case LargeDiag_AWPM:
481: PetscViewerASCIIPrintf(viewer," Row permutation LargeDiag_AWPM\n");
482: break;
483: case MY_PERMR:
484: PetscViewerASCIIPrintf(viewer," Row permutation MY_PERMR\n");
485: break;
486: default:
487: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
488: }
490: switch (options.ColPerm) {
491: case NATURAL:
492: PetscViewerASCIIPrintf(viewer," Column permutation NATURAL\n");
493: break;
494: case MMD_AT_PLUS_A:
495: PetscViewerASCIIPrintf(viewer," Column permutation MMD_AT_PLUS_A\n");
496: break;
497: case MMD_ATA:
498: PetscViewerASCIIPrintf(viewer," Column permutation MMD_ATA\n");
499: break;
500: /* Even though this is called METIS, the SuperLU_DIST code sets this by default if PARMETIS is defined, not METIS */
501: case METIS_AT_PLUS_A:
502: PetscViewerASCIIPrintf(viewer," Column permutation METIS_AT_PLUS_A\n");
503: break;
504: case PARMETIS:
505: PetscViewerASCIIPrintf(viewer," Column permutation PARMETIS\n");
506: break;
507: default:
508: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
509: }
511: PetscViewerASCIIPrintf(viewer," Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);
513: if (lu->FactPattern == SamePattern) {
514: PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern\n");
515: } else if (lu->FactPattern == SamePattern_SameRowPerm) {
516: PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern_SameRowPerm\n");
517: } else if (lu->FactPattern == DOFACT) {
518: PetscViewerASCIIPrintf(viewer," Repeated factorization DOFACT\n");
519: } else {
520: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown factorization pattern");
521: }
522: return(0);
523: }
525: static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
526: {
527: PetscErrorCode ierr;
528: PetscBool iascii;
529: PetscViewerFormat format;
532: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
533: if (iascii) {
534: PetscViewerGetFormat(viewer,&format);
535: if (format == PETSC_VIEWER_ASCII_INFO) {
536: MatView_Info_SuperLU_DIST(A,viewer);
537: }
538: }
539: return(0);
540: }
542: static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
543: {
544: Mat B;
545: Mat_SuperLU_DIST *lu;
546: PetscErrorCode ierr;
547: PetscInt M=A->rmap->N,N=A->cmap->N,indx;
548: PetscMPIInt size;
549: superlu_dist_options_t options;
550: PetscBool flg;
551: const char *colperm[] = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
552: const char *rowperm[] = {"NOROWPERM","LargeDiag_MC64","LargeDiag_AWPM","MY_PERMR"};
553: const char *factPattern[] = {"SamePattern","SamePattern_SameRowPerm","DOFACT"};
554: PetscBool set;
557: /* Create the factorization matrix */
558: MatCreate(PetscObjectComm((PetscObject)A),&B);
559: MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);
560: PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);
561: MatSetUp(B);
562: B->ops->getinfo = MatGetInfo_External;
563: B->ops->view = MatView_SuperLU_DIST;
564: B->ops->destroy = MatDestroy_SuperLU_DIST;
566: /* Set the default input options:
567: options.Fact = DOFACT;
568: options.Equil = YES;
569: options.ParSymbFact = NO;
570: options.ColPerm = METIS_AT_PLUS_A;
571: options.RowPerm = LargeDiag_MC64;
572: options.ReplaceTinyPivot = YES;
573: options.IterRefine = DOUBLE;
574: options.Trans = NOTRANS;
575: options.SolveInitialized = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
576: options.RefineInitialized = NO;
577: options.PrintStat = YES;
578: options.SymPattern = NO;
579: */
580: set_default_options_dist(&options);
582: B->trivialsymbolic = PETSC_TRUE;
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*/