Actual source code: superlu.c
1: #define PETSCMAT_DLL
3: /* --------------------------------------------------------------------
5: This file implements a subclass of the SeqAIJ matrix class that uses
6: the SuperLU 3.0 sparse solver. You can use this as a starting point for
7: implementing your own subclass of a PETSc matrix class.
9: This demonstrates a way to make an implementation inheritence of a PETSc
10: matrix type. This means constructing a new matrix type (SuperLU) by changing some
11: of the methods of the previous type (SeqAIJ), adding additional data, and possibly
12: additional method. (See any book on object oriented programming).
13: */
15: /*
16: Defines the data structure for the base matrix type (SeqAIJ)
17: */
18: #include ../src/mat/impls/aij/seq/aij.h
20: /*
21: SuperLU include files
22: */
24: #if defined(PETSC_USE_COMPLEX)
25: #include "slu_zdefs.h"
26: #else
27: #include "slu_ddefs.h"
28: #endif
29: #include "slu_util.h"
32: /*
33: This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
34: */
35: typedef struct {
36: SuperMatrix A,L,U,B,X;
37: superlu_options_t options;
38: PetscInt *perm_c; /* column permutation vector */
39: PetscInt *perm_r; /* row permutations from partial pivoting */
40: PetscInt *etree;
41: PetscReal *R, *C;
42: char equed[1];
43: PetscInt lwork;
44: void *work;
45: PetscReal rpg, rcond;
46: mem_usage_t mem_usage;
47: MatStructure flg;
49: /* Flag to clean up (non-global) SuperLU objects during Destroy */
50: PetscTruth CleanUpSuperLU;
51: } Mat_SuperLU;
63: /*
64: Utility function
65: */
68: PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
69: {
70: Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr;
71: PetscErrorCode ierr;
72: superlu_options_t options;
75: /* check if matrix is superlu_dist type */
76: if (A->ops->solve != MatSolve_SuperLU) return(0);
78: options = lu->options;
79: PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
80: PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");
81: PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);
82: PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);
83: PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");
84: PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);
85: PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");
86: PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");
87: PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);
88: PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");
89: PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");
90: PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);
91: return(0);
92: }
94: /*
95: These are the methods provided to REPLACE the corresponding methods of the
96: SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
97: */
100: PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
101: {
102: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data;
103: Mat_SuperLU *lu = (Mat_SuperLU*)(F)->spptr;
105: PetscInt sinfo;
106: SuperLUStat_t stat;
107: PetscReal ferr, berr;
108: NCformat *Ustore;
109: SCformat *Lstore;
110:
112: if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
113: lu->options.Fact = SamePattern;
114: /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
115: Destroy_SuperMatrix_Store(&lu->A);
116: if ( lu->lwork >= 0 ) {
117: Destroy_SuperNode_Matrix(&lu->L);
118: Destroy_CompCol_Matrix(&lu->U);
119: lu->options.Fact = SamePattern;
120: }
121: }
123: /* Create the SuperMatrix for lu->A=A^T:
124: Since SuperLU likes column-oriented matrices,we pass it the transpose,
125: and then solve A^T X = B in MatSolve(). */
126: #if defined(PETSC_USE_COMPLEX)
127: zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
128: SLU_NC,SLU_Z,SLU_GE);
129: #else
130: dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
131: SLU_NC,SLU_D,SLU_GE);
132: #endif
133:
134: /* Initialize the statistics variables. */
135: StatInit(&stat);
137: /* Numerical factorization */
138: lu->B.ncol = 0; /* Indicate not to solve the system */
139: #if defined(PETSC_USE_COMPLEX)
140: zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
141: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
142: &lu->mem_usage, &stat, &sinfo);
143: #else
144: dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
145: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
146: &lu->mem_usage, &stat, &sinfo);
147: #endif
148: if ( !sinfo || sinfo == lu->A.ncol+1 ) {
149: if ( lu->options.PivotGrowth )
150: PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg);
151: if ( lu->options.ConditionNumber )
152: PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond);
153: } else if ( sinfo > 0 ){
154: if ( lu->lwork == -1 ) {
155: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
156: } else {
157: PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo);
158: }
159: } else { /* sinfo < 0 */
160: SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
161: }
163: if ( lu->options.PrintStat ) {
164: PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
165: StatPrint(&stat);
166: Lstore = (SCformat *) lu->L.Store;
167: Ustore = (NCformat *) lu->U.Store;
168: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz);
169: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz);
170: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
171: PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n",
172: lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
173: lu->mem_usage.expansions);
174: }
175: StatFree(&stat);
177: lu->flg = SAME_NONZERO_PATTERN;
178: (F)->ops->solve = MatSolve_SuperLU;
179: (F)->ops->solvetranspose = MatSolveTranspose_SuperLU;
180: return(0);
181: }
185: PetscErrorCode MatDestroy_SuperLU(Mat A)
186: {
188: Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr;
191: if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
192: Destroy_SuperMatrix_Store(&lu->A);
193: Destroy_SuperMatrix_Store(&lu->B);
194: Destroy_SuperMatrix_Store(&lu->X);
196: PetscFree(lu->etree);
197: PetscFree(lu->perm_r);
198: PetscFree(lu->perm_c);
199: PetscFree(lu->R);
200: PetscFree(lu->C);
201: if ( lu->lwork >= 0 ) {
202: Destroy_SuperNode_Matrix(&lu->L);
203: Destroy_CompCol_Matrix(&lu->U);
204: }
205: }
206: MatDestroy_SeqAIJ(A);
207: return(0);
208: }
212: PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
213: {
214: PetscErrorCode ierr;
215: PetscTruth iascii;
216: PetscViewerFormat format;
219: MatView_SeqAIJ(A,viewer);
221: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
222: if (iascii) {
223: PetscViewerGetFormat(viewer,&format);
224: if (format == PETSC_VIEWER_ASCII_INFO) {
225: MatFactorInfo_SuperLU(A,viewer);
226: }
227: }
228: return(0);
229: }
234: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
235: {
236: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
237: PetscScalar *barray,*xarray;
239: PetscInt info,i;
240: SuperLUStat_t stat;
241: PetscReal ferr,berr;
244: if ( lu->lwork == -1 ) {
245: return(0);
246: }
247: lu->B.ncol = 1; /* Set the number of right-hand side */
248: VecGetArray(b,&barray);
249: VecGetArray(x,&xarray);
251: #if defined(PETSC_USE_COMPLEX)
252: ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
253: ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
254: #else
255: ((DNformat*)lu->B.Store)->nzval = barray;
256: ((DNformat*)lu->X.Store)->nzval = xarray;
257: #endif
259: /* Initialize the statistics variables. */
260: StatInit(&stat);
262: lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
263: #if defined(PETSC_USE_COMPLEX)
264: zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
265: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
266: &lu->mem_usage, &stat, &info);
267: #else
268: dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
269: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
270: &lu->mem_usage, &stat, &info);
271: #endif
272: VecRestoreArray(b,&barray);
273: VecRestoreArray(x,&xarray);
275: if ( !info || info == lu->A.ncol+1 ) {
276: if ( lu->options.IterRefine ) {
277: PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
278: PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
279: for (i = 0; i < 1; ++i)
280: PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
281: }
282: } else if ( info > 0 ){
283: if ( lu->lwork == -1 ) {
284: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol);
285: } else {
286: PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info);
287: }
288: } else if (info < 0){
289: SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
290: }
292: if ( lu->options.PrintStat ) {
293: PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
294: StatPrint(&stat);
295: }
296: StatFree(&stat);
297: return(0);
298: }
302: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
303: {
304: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
308: lu->options.Trans = TRANS;
309: MatSolve_SuperLU_Private(A,b,x);
310: return(0);
311: }
315: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
316: {
317: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
321: lu->options.Trans = NOTRANS;
322: MatSolve_SuperLU_Private(A,b,x);
323: return(0);
324: }
327: /*
328: Note the r permutation is ignored
329: */
332: PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
333: {
334: Mat_SuperLU *lu = (Mat_SuperLU*)((F)->spptr);
336: PetscInt m=A->rmap->n,n=A->cmap->n;
340: /* Allocate spaces (notice sizes are for the transpose) */
341: PetscMalloc(m*sizeof(PetscInt),&lu->etree);
342: PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);
343: PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);
344: PetscMalloc(n*sizeof(PetscInt),&lu->R);
345: PetscMalloc(m*sizeof(PetscInt),&lu->C);
346:
347: /* create rhs and solution x without allocate space for .Store */
348: #if defined(PETSC_USE_COMPLEX)
349: zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
350: zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
351: #else
352: dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
353: dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
354: #endif
356: lu->flg = DIFFERENT_NONZERO_PATTERN;
357: lu->CleanUpSuperLU = PETSC_TRUE;
358: (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
359: return(0);
360: }
365: PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
366: {
368: *type = MAT_SOLVER_SUPERLU;
369: return(0);
370: }
372:
374: /*MC
375: MAT_SOLVER_SUPERLU = "superlu" - A solver package roviding direct solvers (LU) for sequential matrices
376: via the external package SuperLU.
378: Use config/configure.py --download-superlu to have PETSc installed with SuperLU
380: Options Database Keys:
381: + -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
382: 1: MMD applied to A'*A,
383: 2: MMD applied to A'+A,
384: 3: COLAMD, approximate minimum degree column ordering
385: . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
386: choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
387: - -mat_superlu_printstat - print SuperLU statistics about the factorization
389: Notes: Do not confuse this with MAT_SOLVER_SUPERLU_DIST which is for parallel sparse solves
391: Level: beginner
393: .seealso: PCLU, MAT_SOLVER_SUPERLU_DIST, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
394: M*/
399: PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
400: {
401: Mat B;
402: Mat_SuperLU *lu;
404: PetscInt indx;
405: PetscTruth flg;
406: const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
407: const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
408: const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
411: MatCreate(((PetscObject)A)->comm,&B);
412: MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
413: MatSetType(B,((PetscObject)A)->type_name);
414: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
416: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
417: B->ops->destroy = MatDestroy_SuperLU;
418: B->ops->view = MatView_SuperLU;
419: B->factor = MAT_FACTOR_LU;
420: B->assembled = PETSC_TRUE; /* required by -ksp_view */
421: B->preallocated = PETSC_TRUE;
422:
423: PetscNewLog(B,Mat_SuperLU,&lu);
424: set_default_options(&lu->options);
425: /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
426: lu->options.Equil = NO;
427: lu->options.PrintStat = NO;
428: lu->lwork = 0; /* allocate space internally by system malloc */
430: PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");
431: PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
432: if (flg) {lu->options.ColPerm = (colperm_t)indx;}
433: PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
434: if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
435: PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);
436: if (flg) lu->options.SymmetricMode = YES;
437: PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);
438: PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);
439: if (flg) lu->options.PivotGrowth = YES;
440: PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);
441: if (flg) lu->options.ConditionNumber = YES;
442: PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);
443: if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
444: PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);
445: if (flg) lu->options.ReplaceTinyPivot = YES;
446: PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);
447: if (flg) lu->options.PrintStat = YES;
448: PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);
449: if (lu->lwork > 0 ){
450: PetscMalloc(lu->lwork,&lu->work);
451: } else if (lu->lwork != 0 && lu->lwork != -1){
452: PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
453: lu->lwork = 0;
454: }
455: PetscOptionsEnd();
457: #ifdef SUPERLU2
458: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);
459: #endif
460: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);
461: B->spptr = lu;
462: *F = B;
463: return(0);
464: }