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 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;
48: SuperLUStat_t stat;
50: /* Flag to clean up (non-global) SuperLU objects during Destroy */
51: PetscTruth CleanUpSuperLU;
52: } Mat_SuperLU;
64: /*
65: Utility function
66: */
69: PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
70: {
71: Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr;
72: PetscErrorCode ierr;
73: superlu_options_t options;
76: /* check if matrix is superlu_dist type */
77: if (A->ops->solve != MatSolve_SuperLU) return(0);
79: options = lu->options;
80: PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
81: PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");
82: PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);
83: PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);
84: PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");
85: PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);
86: PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");
87: PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");
88: PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);
89: PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");
90: PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");
91: PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);
92: if (A->factor == MAT_FACTOR_ILU){
93: PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);
94: PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);
95: PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);
96: PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);
97: PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);
98: PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);
99: }
100: return(0);
101: }
103: /*
104: These are the methods provided to REPLACE the corresponding methods of the
105: SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
106: */
109: PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
110: {
111: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data;
112: Mat_SuperLU *lu = (Mat_SuperLU*)(F)->spptr;
114: PetscInt sinfo;
115: PetscReal ferr, berr;
116: NCformat *Ustore;
117: SCformat *Lstore;
118:
120: if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
121: lu->options.Fact = SamePattern;
122: /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
123: Destroy_SuperMatrix_Store(&lu->A);
124: if ( lu->lwork >= 0 ) {
125: Destroy_SuperNode_Matrix(&lu->L);
126: Destroy_CompCol_Matrix(&lu->U);
127: lu->options.Fact = SamePattern;
128: }
129: }
131: /* Create the SuperMatrix for lu->A=A^T:
132: Since SuperLU likes column-oriented matrices,we pass it the transpose,
133: and then solve A^T X = B in MatSolve(). */
134: #if defined(PETSC_USE_COMPLEX)
135: zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
136: SLU_NC,SLU_Z,SLU_GE);
137: #else
138: dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
139: SLU_NC,SLU_D,SLU_GE);
140: #endif
142: /* Numerical factorization */
143: lu->B.ncol = 0; /* Indicate not to solve the system */
144: if (F->factor == MAT_FACTOR_LU){
145: #if defined(PETSC_USE_COMPLEX)
146: zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
147: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
148: &lu->mem_usage, &lu->stat, &sinfo);
149: #else
150: dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
151: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
152: &lu->mem_usage, &lu->stat, &sinfo);
153: #endif
154: } else if (F->factor == MAT_FACTOR_ILU){
155: /* Compute the incomplete factorization, condition number and pivot growth */
156: #if defined(PETSC_USE_COMPLEX)
157: zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
158: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
159: &lu->mem_usage, &lu->stat, &sinfo);
160: #else
161: dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
162: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
163: &lu->mem_usage, &lu->stat, &sinfo);
164: #endif
165: } else {
166: SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
167: }
168: if ( !sinfo || sinfo == lu->A.ncol+1 ) {
169: if ( lu->options.PivotGrowth )
170: PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg);
171: if ( lu->options.ConditionNumber )
172: PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond);
173: } else if ( sinfo > 0 ){
174: if ( lu->lwork == -1 ) {
175: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
176: } else {
177: PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo);
178: }
179: } else { /* sinfo < 0 */
180: SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
181: }
183: if ( lu->options.PrintStat ) {
184: PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
185: StatPrint(&lu->stat);
186: Lstore = (SCformat *) lu->L.Store;
187: Ustore = (NCformat *) lu->U.Store;
188: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz);
189: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz);
190: PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
191: PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n",
192: lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
193: }
195: lu->flg = SAME_NONZERO_PATTERN;
196: (F)->ops->solve = MatSolve_SuperLU;
197: (F)->ops->solvetranspose = MatSolveTranspose_SuperLU;
198: return(0);
199: }
203: PetscErrorCode MatDestroy_SuperLU(Mat A)
204: {
206: Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr;
209: if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
210: Destroy_SuperMatrix_Store(&lu->A);
211: Destroy_SuperMatrix_Store(&lu->B);
212: Destroy_SuperMatrix_Store(&lu->X);
213: StatFree(&lu->stat);
215: PetscFree(lu->etree);
216: PetscFree(lu->perm_r);
217: PetscFree(lu->perm_c);
218: PetscFree(lu->R);
219: PetscFree(lu->C);
220: if ( lu->lwork >= 0 ) {
221: Destroy_SuperNode_Matrix(&lu->L);
222: Destroy_CompCol_Matrix(&lu->U);
223: }
224: }
225: MatDestroy_SeqAIJ(A);
226: return(0);
227: }
231: PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
232: {
233: PetscErrorCode ierr;
234: PetscTruth iascii;
235: PetscViewerFormat format;
238: MatView_SeqAIJ(A,viewer);
240: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
241: if (iascii) {
242: PetscViewerGetFormat(viewer,&format);
243: if (format == PETSC_VIEWER_ASCII_INFO) {
244: MatFactorInfo_SuperLU(A,viewer);
245: }
246: }
247: return(0);
248: }
253: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
254: {
255: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
256: PetscScalar *barray,*xarray;
258: PetscInt info,i;
259: PetscReal ferr,berr;
262: if ( lu->lwork == -1 ) {
263: return(0);
264: }
265: lu->B.ncol = 1; /* Set the number of right-hand side */
266: VecGetArray(b,&barray);
267: VecGetArray(x,&xarray);
269: #if defined(PETSC_USE_COMPLEX)
270: ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
271: ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
272: #else
273: ((DNformat*)lu->B.Store)->nzval = barray;
274: ((DNformat*)lu->X.Store)->nzval = xarray;
275: #endif
277: lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
278: if (A->factor == MAT_FACTOR_LU){
279: #if defined(PETSC_USE_COMPLEX)
280: zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
281: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
282: &lu->mem_usage, &lu->stat, &info);
283: #else
284: dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
285: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
286: &lu->mem_usage, &lu->stat, &info);
287: #endif
288: } else if (A->factor == MAT_FACTOR_ILU){
289: #if defined(PETSC_USE_COMPLEX)
290: zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
291: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
292: &lu->mem_usage, &lu->stat, &info);
293: #else
294: dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
295: &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
296: &lu->mem_usage, &lu->stat, &info);
297: #endif
298: } else {
299: SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
300: }
301: VecRestoreArray(b,&barray);
302: VecRestoreArray(x,&xarray);
304: if ( !info || info == lu->A.ncol+1 ) {
305: if ( lu->options.IterRefine ) {
306: PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
307: PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
308: for (i = 0; i < 1; ++i)
309: PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
310: }
311: } else if ( info > 0 ){
312: if ( lu->lwork == -1 ) {
313: PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol);
314: } else {
315: PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info);
316: }
317: } else if (info < 0){
318: SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
319: }
321: if ( lu->options.PrintStat ) {
322: PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
323: StatPrint(&lu->stat);
324: }
325: return(0);
326: }
330: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
331: {
332: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
336: lu->options.Trans = TRANS;
337: MatSolve_SuperLU_Private(A,b,x);
338: return(0);
339: }
343: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
344: {
345: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
349: lu->options.Trans = NOTRANS;
350: MatSolve_SuperLU_Private(A,b,x);
351: return(0);
352: }
354: /*
355: Note the r permutation is ignored
356: */
359: PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
360: {
361: Mat_SuperLU *lu = (Mat_SuperLU*)((F)->spptr);
362:
364: lu->flg = DIFFERENT_NONZERO_PATTERN;
365: lu->CleanUpSuperLU = PETSC_TRUE;
366: (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
367: return(0);
368: }
373: PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
374: {
376: *type = MAT_SOLVER_SUPERLU;
377: return(0);
378: }
380:
382: /*MC
383: MAT_SOLVER_SUPERLU = "superlu" - A solver package roviding direct solvers (LU) for sequential matrices
384: via the external package SuperLU.
386: Use config/configure.py --download-superlu to have PETSc installed with SuperLU
388: Options Database Keys:
389: + -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
390: 1: MMD applied to A'*A,
391: 2: MMD applied to A'+A,
392: 3: COLAMD, approximate minimum degree column ordering
393: . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
394: choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
395: - -mat_superlu_printstat - print SuperLU statistics about the factorization
397: Notes: Do not confuse this with MAT_SOLVER_SUPERLU_DIST which is for parallel sparse solves
399: Level: beginner
401: .seealso: PCLU, MAT_SOLVER_SUPERLU_DIST, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
402: M*/
407: PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
408: {
409: Mat B;
410: Mat_SuperLU *lu;
412: PetscInt indx,m=A->rmap->n,n=A->cmap->n;
413: PetscTruth flg;
414: const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
415: const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
416: const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
419: MatCreate(((PetscObject)A)->comm,&B);
420: MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
421: MatSetType(B,((PetscObject)A)->type_name);
422: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
424: if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
425: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
426: B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
427: } else {
428: SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
429: }
431: B->ops->destroy = MatDestroy_SuperLU;
432: B->ops->view = MatView_SuperLU;
433: B->factor = ftype;
434: B->assembled = PETSC_TRUE; /* required by -ksp_view */
435: B->preallocated = PETSC_TRUE;
436:
437: PetscNewLog(B,Mat_SuperLU,&lu);
438: if (ftype == MAT_FACTOR_LU){
439: set_default_options(&lu->options);
440: } else if (ftype == MAT_FACTOR_ILU){
441: /* Set the default input options of ilu:
442: options.Fact = DOFACT;
443: options.Equil = YES;
444: options.ColPerm = COLAMD;
445: options.DiagPivotThresh = 0.1; //different from complete LU
446: options.Trans = NOTRANS;
447: options.IterRefine = NOREFINE;
448: options.SymmetricMode = NO;
449: options.PivotGrowth = NO;
450: options.ConditionNumber = NO;
451: options.PrintStat = YES;
452: options.RowPerm = LargeDiag;
453: options.ILU_DropTol = 1e-4;
454: options.ILU_FillTol = 1e-2;
455: options.ILU_FillFactor = 10.0;
456: options.ILU_DropRule = DROP_BASIC | DROP_AREA;
457: options.ILU_Norm = INF_NORM;
458: options.ILU_MILU = SMILU_2;
459: */
460:
461: ilu_set_default_options(&lu->options);
462: }
463: /* equilibration causes error in solve()(ref. [petsc-maint #42782] SuperLU troubles)
464: thus not supported here. See dgssvx.c for possible reason. */
465: lu->options.Equil = NO;
466: lu->options.PrintStat = NO;
468: /* Initialize the statistics variables. */
469: StatInit(&lu->stat);
470: lu->lwork = 0; /* allocate space internally by system malloc */
472: PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");
473: PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);
474: if (flg) lu->options.Equil = YES;
475: PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
476: if (flg) {lu->options.ColPerm = (colperm_t)indx;}
477: PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
478: if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
479: PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);
480: if (flg) lu->options.SymmetricMode = YES;
481: PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);
482: PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);
483: if (flg) lu->options.PivotGrowth = YES;
484: PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);
485: if (flg) lu->options.ConditionNumber = YES;
486: PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);
487: if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
488: PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);
489: if (flg) lu->options.ReplaceTinyPivot = YES;
490: PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);
491: if (flg) lu->options.PrintStat = YES;
492: PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);
493: if (lu->lwork > 0 ){
494: PetscMalloc(lu->lwork,&lu->work);
495: } else if (lu->lwork != 0 && lu->lwork != -1){
496: PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
497: lu->lwork = 0;
498: }
499: /* ilu options */
500: PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);
501: PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);
502: PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);
503: PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);
504: PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);
505: if (flg){
506: lu->options.ILU_Norm = (norm_t)indx;
507: }
508: PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);
509: if (flg){
510: lu->options.ILU_MILU = (milu_t)indx;
511: }
512: PetscOptionsEnd();
514: /* Allocate spaces (notice sizes are for the transpose) */
515: PetscMalloc(m*sizeof(PetscInt),&lu->etree);
516: PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);
517: PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);
518: PetscMalloc(n*sizeof(PetscScalar),&lu->R);
519: PetscMalloc(m*sizeof(PetscScalar),&lu->C);
520:
521: /* create rhs and solution x without allocate space for .Store */
522: #if defined(PETSC_USE_COMPLEX)
523: zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
524: zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
525: #else
526: dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
527: dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
528: #endif
530: #ifdef SUPERLU2
531: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);
532: #endif
533: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);
534: B->spptr = lu;
535: *F = B;
536: return(0);
537: }