Actual source code: superlu.c
petsc-3.6.1 2015-08-06
2: /* --------------------------------------------------------------------
4: This file implements a subclass of the SeqAIJ matrix class that uses
5: the SuperLU sparse solver.
6: */
8: /*
9: Defines the data structure for the base matrix type (SeqAIJ)
10: */
11: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
13: /*
14: SuperLU include files
15: */
16: EXTERN_C_BEGIN
17: #if defined(PETSC_USE_COMPLEX)
18: #if defined(PETSC_USE_REAL_SINGLE)
19: #include <slu_cdefs.h>
20: #else
21: #include <slu_zdefs.h>
22: #endif
23: #else
24: #if defined(PETSC_USE_REAL_SINGLE)
25: #include <slu_sdefs.h>
26: #else
27: #include <slu_ddefs.h>
28: #endif
29: #endif
30: #include <slu_util.h>
31: EXTERN_C_END
33: /*
34: This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
35: */
36: typedef struct {
37: SuperMatrix A,L,U,B,X;
38: superlu_options_t options;
39: PetscInt *perm_c; /* column permutation vector */
40: PetscInt *perm_r; /* row permutations from partial pivoting */
41: PetscInt *etree;
42: PetscReal *R, *C;
43: char equed[1];
44: PetscInt lwork;
45: void *work;
46: PetscReal rpg, rcond;
47: mem_usage_t mem_usage;
48: MatStructure flg;
49: SuperLUStat_t stat;
50: Mat A_dup;
51: PetscScalar *rhs_dup;
53: /* Flag to clean up (non-global) SuperLU objects during Destroy */
54: PetscBool CleanUpSuperLU;
55: } Mat_SuperLU;
57: extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
58: extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo*);
59: extern PetscErrorCode MatDestroy_SuperLU(Mat);
60: extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
61: extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
62: extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
63: extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat);
64: extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
65: extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*);
66: extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat*);
68: /*
69: Utility function
70: */
73: PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
74: {
75: Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr;
76: PetscErrorCode ierr;
77: superlu_options_t options;
80: /* check if matrix is superlu_dist type */
81: if (A->ops->solve != MatSolve_SuperLU) return(0);
83: options = lu->options;
85: PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
86: PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");
87: PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);
88: PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);
89: PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");
90: PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);
91: PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");
92: PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");
93: PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);
94: PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");
95: PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");
96: PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);
97: if (A->factortype == MAT_FACTOR_ILU) {
98: PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);
99: PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);
100: PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);
101: PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);
102: PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);
103: PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);
104: }
105: return(0);
106: }
108: /*
109: These are the methods provided to REPLACE the corresponding methods of the
110: SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
111: */
114: PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
115: {
116: Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr;
117: Mat_SeqAIJ *aa;
119: PetscInt sinfo;
120: PetscReal ferr, berr;
121: NCformat *Ustore;
122: SCformat *Lstore;
125: if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */
126: lu->options.Fact = SamePattern;
127: /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
128: Destroy_SuperMatrix_Store(&lu->A);
129: if (lu->options.Equil) {
130: MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);
131: }
132: if (lu->lwork >= 0) {
133: PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
134: PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
135: lu->options.Fact = SamePattern;
136: }
137: }
139: /* Create the SuperMatrix for lu->A=A^T:
140: Since SuperLU likes column-oriented matrices,we pass it the transpose,
141: and then solve A^T X = B in MatSolve(). */
142: if (lu->options.Equil) {
143: aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
144: } else {
145: aa = (Mat_SeqAIJ*)(A)->data;
146: }
147: #if defined(PETSC_USE_COMPLEX)
148: #if defined(PETSC_USE_REAL_SINGLE)
149: PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE));
150: #else
151: PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE));
152: #endif
153: #else
154: #if defined(PETSC_USE_REAL_SINGLE)
155: PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE));
156: #else
157: PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE));
158: #endif
159: #endif
161: /* Numerical factorization */
162: lu->B.ncol = 0; /* Indicate not to solve the system */
163: if (F->factortype == MAT_FACTOR_LU) {
164: #if defined(PETSC_USE_COMPLEX)
165: #if defined(PETSC_USE_REAL_SINGLE)
166: PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
167: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
168: &lu->mem_usage, &lu->stat, &sinfo));
169: #else
170: PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
171: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
172: &lu->mem_usage, &lu->stat, &sinfo));
173: #endif
174: #else
175: #if defined(PETSC_USE_REAL_SINGLE)
176: PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
177: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
178: &lu->mem_usage, &lu->stat, &sinfo));
179: #else
180: PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
181: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
182: &lu->mem_usage, &lu->stat, &sinfo));
183: #endif
184: #endif
185: } else if (F->factortype == MAT_FACTOR_ILU) {
186: /* Compute the incomplete factorization, condition number and pivot growth */
187: #if defined(PETSC_USE_COMPLEX)
188: #if defined(PETSC_USE_REAL_SINGLE)
189: PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
190: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
191: &lu->mem_usage, &lu->stat, &sinfo));
192: #else
193: PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
194: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
195: &lu->mem_usage, &lu->stat, &sinfo));
196: #endif
197: #else
198: #if defined(PETSC_USE_REAL_SINGLE)
199: PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
200: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
201: &lu->mem_usage, &lu->stat, &sinfo));
202: #else
203: PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
204: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
205: &lu->mem_usage, &lu->stat, &sinfo));
206: #endif
207: #endif
208: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
209: if (!sinfo || sinfo == lu->A.ncol+1) {
210: if (lu->options.PivotGrowth) {
211: PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg);
212: }
213: if (lu->options.ConditionNumber) {
214: PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond);
215: }
216: } else if (sinfo > 0) {
217: if (lu->lwork == -1) {
218: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
219: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
220: } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
222: if (lu->options.PrintStat) {
223: PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
224: PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
225: Lstore = (SCformat*) lu->L.Store;
226: Ustore = (NCformat*) lu->U.Store;
227: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz);
228: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz);
229: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
230: PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n",
231: lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
232: }
234: lu->flg = SAME_NONZERO_PATTERN;
235: F->ops->solve = MatSolve_SuperLU;
236: F->ops->solvetranspose = MatSolveTranspose_SuperLU;
237: F->ops->matsolve = NULL;
238: return(0);
239: }
243: PetscErrorCode MatGetDiagonal_SuperLU(Mat A,Vec v)
244: {
246: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: SuperLU factor");
247: return(0);
248: }
252: PetscErrorCode MatDestroy_SuperLU(Mat A)
253: {
255: Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr;
258: if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
259: PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A));
260: PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B));
261: PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X));
262: PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat));
263: if (lu->lwork >= 0) {
264: PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
265: PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
266: }
267: }
268: if (lu) {
269: PetscFree(lu->etree);
270: PetscFree(lu->perm_r);
271: PetscFree(lu->perm_c);
272: PetscFree(lu->R);
273: PetscFree(lu->C);
274: PetscFree(lu->rhs_dup);
275: MatDestroy(&lu->A_dup);
276: }
277: PetscFree(A->spptr);
279: /* clear composed functions */
280: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
281: PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);
283: MatDestroy_SeqAIJ(A);
284: return(0);
285: }
289: PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
290: {
291: PetscErrorCode ierr;
292: PetscBool iascii;
293: PetscViewerFormat format;
296: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
297: if (iascii) {
298: PetscViewerGetFormat(viewer,&format);
299: if (format == PETSC_VIEWER_ASCII_INFO) {
300: MatFactorInfo_SuperLU(A,viewer);
301: }
302: }
303: return(0);
304: }
309: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
310: {
311: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
312: const PetscScalar *barray;
313: PetscScalar *xarray;
314: PetscErrorCode ierr;
315: PetscInt info,i,n;
316: PetscReal ferr,berr;
317: static PetscBool cite = PETSC_FALSE;
320: if (lu->lwork == -1) return(0);
321: PetscCitationsRegister("@article{superlu99,\n author = {James W. Demmel and Stanley C. Eisenstat and\n John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n title = {A supernodal approach to sparse partial pivoting},\n journal = {SIAM J. Matrix Analysis and Applications},\n year = {1999},\n volume = {20},\n number = {3},\n pages = {720-755}\n}\n",&cite);
323: VecGetLocalSize(x,&n);
324: lu->B.ncol = 1; /* Set the number of right-hand side */
325: if (lu->options.Equil && !lu->rhs_dup) {
326: /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
327: PetscMalloc1(n,&lu->rhs_dup);
328: }
329: if (lu->options.Equil) {
330: /* Copy b into rsh_dup */
331: VecGetArrayRead(b,&barray);
332: PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));
333: VecRestoreArrayRead(b,&barray);
334: barray = lu->rhs_dup;
335: } else {
336: VecGetArrayRead(b,&barray);
337: }
338: VecGetArray(x,&xarray);
340: #if defined(PETSC_USE_COMPLEX)
341: #if defined(PETSC_USE_REAL_SINGLE)
342: ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
343: ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
344: #else
345: ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
346: ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
347: #endif
348: #else
349: ((DNformat*)lu->B.Store)->nzval = (void*)barray;
350: ((DNformat*)lu->X.Store)->nzval = xarray;
351: #endif
353: lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
354: if (A->factortype == MAT_FACTOR_LU) {
355: #if defined(PETSC_USE_COMPLEX)
356: #if defined(PETSC_USE_REAL_SINGLE)
357: PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
358: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
359: &lu->mem_usage, &lu->stat, &info));
360: #else
361: PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
362: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
363: &lu->mem_usage, &lu->stat, &info));
364: #endif
365: #else
366: #if defined(PETSC_USE_REAL_SINGLE)
367: PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
368: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
369: &lu->mem_usage, &lu->stat, &info));
370: #else
371: PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
372: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
373: &lu->mem_usage, &lu->stat, &info));
374: #endif
375: #endif
376: } else if (A->factortype == MAT_FACTOR_ILU) {
377: #if defined(PETSC_USE_COMPLEX)
378: #if defined(PETSC_USE_REAL_SINGLE)
379: PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
380: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
381: &lu->mem_usage, &lu->stat, &info));
382: #else
383: PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
384: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
385: &lu->mem_usage, &lu->stat, &info));
386: #endif
387: #else
388: #if defined(PETSC_USE_REAL_SINGLE)
389: PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
390: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
391: &lu->mem_usage, &lu->stat, &info));
392: #else
393: PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
394: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
395: &lu->mem_usage, &lu->stat, &info));
396: #endif
397: #endif
398: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
399: if (!lu->options.Equil) {
400: VecRestoreArrayRead(b,&barray);
401: }
402: VecRestoreArray(x,&xarray);
404: if (!info || info == lu->A.ncol+1) {
405: if (lu->options.IterRefine) {
406: PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
407: PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
408: for (i = 0; i < 1; ++i) {
409: PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
410: }
411: }
412: } else if (info > 0) {
413: if (lu->lwork == -1) {
414: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol);
415: } else {
416: PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info);
417: }
418: } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
420: if (lu->options.PrintStat) {
421: PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
422: PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
423: }
424: return(0);
425: }
429: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
430: {
431: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
435: lu->options.Trans = TRANS;
437: MatSolve_SuperLU_Private(A,b,x);
438: return(0);
439: }
443: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
444: {
445: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
449: lu->options.Trans = NOTRANS;
451: MatSolve_SuperLU_Private(A,b,x);
452: return(0);
453: }
457: PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
458: {
459: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
460: PetscBool flg;
464: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
465: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
466: PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
467: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
468: lu->options.Trans = TRANS;
469: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
470: return(0);
471: }
473: /*
474: Note the r permutation is ignored
475: */
478: PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
479: {
480: Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr);
483: lu->flg = DIFFERENT_NONZERO_PATTERN;
484: lu->CleanUpSuperLU = PETSC_TRUE;
485: F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
486: return(0);
487: }
491: static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
492: {
493: Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
496: lu->options.ILU_DropTol = dtol;
497: return(0);
498: }
502: /*@
503: MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
504: Logically Collective on Mat
506: Input Parameters:
507: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
508: - dtol - drop tolerance
510: Options Database:
511: . -mat_superlu_ilu_droptol <dtol>
513: Level: beginner
515: References: SuperLU Users' Guide
517: .seealso: MatGetFactor()
518: @*/
519: PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
520: {
526: PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));
527: return(0);
528: }
532: PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
533: {
535: *type = MATSOLVERSUPERLU;
536: return(0);
537: }
540: /*MC
541: MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
542: via the external package SuperLU.
544: Use ./configure --download-superlu to have PETSc installed with SuperLU
546: Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver
548: Options Database Keys:
549: + -mat_superlu_equil <FALSE> - Equil (None)
550: . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
551: . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
552: . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None)
553: . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None)
554: . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None)
555: . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None)
556: . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag
557: . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
558: . -mat_superlu_printstat <FALSE> - PrintStat (None)
559: . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None)
560: . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None)
561: . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None)
562: . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None)
563: . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None)
564: . -mat_superlu_ilu_norm <0> - ILU_Norm (None)
565: - -mat_superlu_ilu_milu <0> - ILU_MILU (None)
567: Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
569: Level: beginner
571: .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
572: M*/
576: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
577: {
578: Mat B;
579: Mat_SuperLU *lu;
581: PetscInt indx,m=A->rmap->n,n=A->cmap->n;
582: PetscBool flg,set;
583: PetscReal real_input;
584: const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
585: const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
586: const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
589: MatCreate(PetscObjectComm((PetscObject)A),&B);
590: MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
591: MatSetType(B,((PetscObject)A)->type_name);
592: MatSeqAIJSetPreallocation(B,0,NULL);
594: if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
595: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
596: B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
597: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
599: B->ops->destroy = MatDestroy_SuperLU;
600: B->ops->view = MatView_SuperLU;
601: B->ops->getdiagonal = MatGetDiagonal_SuperLU;
602: B->factortype = ftype;
603: B->assembled = PETSC_TRUE; /* required by -ksp_view */
604: B->preallocated = PETSC_TRUE;
606: PetscNewLog(B,&lu);
608: if (ftype == MAT_FACTOR_LU) {
609: set_default_options(&lu->options);
610: /* Comments from SuperLU_4.0/SRC/dgssvx.c:
611: "Whether or not the system will be equilibrated depends on the
612: scaling of the matrix A, but if equilibration is used, A is
613: overwritten by diag(R)*A*diag(C) and B by diag(R)*B
614: (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
615: We set 'options.Equil = NO' as default because additional space is needed for it.
616: */
617: lu->options.Equil = NO;
618: } else if (ftype == MAT_FACTOR_ILU) {
619: /* Set the default input options of ilu: */
620: PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options));
621: }
622: lu->options.PrintStat = NO;
624: /* Initialize the statistics variables. */
625: PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat));
626: lu->lwork = 0; /* allocate space internally by system malloc */
628: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");
629: PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);
630: PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
631: if (flg) lu->options.ColPerm = (colperm_t)indx;
632: PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
633: if (flg) lu->options.IterRefine = (IterRefine_t)indx;
634: PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);
635: if (set && flg) lu->options.SymmetricMode = YES;
636: PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);
637: if (flg) lu->options.DiagPivotThresh = (double) real_input;
638: PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);
639: if (set && flg) lu->options.PivotGrowth = YES;
640: PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);
641: if (set && flg) lu->options.ConditionNumber = YES;
642: PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);
643: if (flg) lu->options.RowPerm = (rowperm_t)indx;
644: PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);
645: if (set && flg) lu->options.ReplaceTinyPivot = YES;
646: PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);
647: if (set && flg) lu->options.PrintStat = YES;
648: PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);
649: if (lu->lwork > 0) {
650: /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/
651: PetscMalloc(lu->lwork,&lu->work);
652: } else if (lu->lwork != 0 && lu->lwork != -1) {
653: PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
654: lu->lwork = 0;
655: }
656: /* ilu options */
657: PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);
658: if (flg) lu->options.ILU_DropTol = (double) real_input;
659: PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);
660: if (flg) lu->options.ILU_FillTol = (double) real_input;
661: PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);
662: if (flg) lu->options.ILU_FillFactor = (double) real_input;
663: PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);
664: PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);
665: if (flg) lu->options.ILU_Norm = (norm_t)indx;
666: PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);
667: if (flg) lu->options.ILU_MILU = (milu_t)indx;
668: PetscOptionsEnd();
669: if (lu->options.Equil == YES) {
670: /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
671: MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);
672: }
674: /* Allocate spaces (notice sizes are for the transpose) */
675: PetscMalloc1(m,&lu->etree);
676: PetscMalloc1(n,&lu->perm_r);
677: PetscMalloc1(m,&lu->perm_c);
678: PetscMalloc1(n,&lu->R);
679: PetscMalloc1(m,&lu->C);
681: /* create rhs and solution x without allocate space for .Store */
682: #if defined(PETSC_USE_COMPLEX)
683: #if defined(PETSC_USE_REAL_SINGLE)
684: PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
685: PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
686: #else
687: PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
688: PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
689: #endif
690: #else
691: #if defined(PETSC_USE_REAL_SINGLE)
692: PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
693: PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
694: #else
695: PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
696: PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
697: #endif
698: #endif
700: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);
701: PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);
702: B->spptr = lu;
704: *F = B;
705: return(0);
706: }
710: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void)
711: {
715: MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);
716: MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);
717: return(0);
718: }