xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision a4efd8ea4bb0600a5aef92b2c23cbfa314827193)
1 #define PETSCMAT_DLL
2 
3 /*  --------------------------------------------------------------------
4 
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.
8 
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 */
14 
15 /*
16      Defines the data structure for the base matrix type (SeqAIJ)
17 */
18 #include "src/mat/impls/aij/seq/aij.h"
19 
20 /*
21      SuperLU include files
22 */
23 EXTERN_C_BEGIN
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"
30 EXTERN_C_END
31 
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 
49   /* Flag to clean up (non-global) SuperLU objects during Destroy */
50   PetscTruth CleanUpSuperLU;
51 } Mat_SuperLU;
52 
53 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
54 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,MatFactorInfo *,Mat *);
55 extern PetscErrorCode MatDestroy_SuperLU(Mat);
56 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
57 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
58 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
59 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
60 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo *,Mat *);
61 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *);
62 
63 /*
64     Utility function
65 */
66 #undef __FUNCT__
67 #define __FUNCT__ "MatFactorInfo_SuperLU"
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;
73 
74   PetscFunctionBegin;
75   /* check if matrix is superlu_dist type */
76   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
77 
78   options = lu->options;
79   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
80   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
81   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
82   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
83   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
84   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
85   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
86   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
87   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
88   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
89   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
90   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
91 
92   PetscFunctionReturn(0);
93 }
94 
95 /*
96     These are the methods provided to REPLACE the corresponding methods of the
97    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
98 */
99 #undef __FUNCT__
100 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
101 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F)
102 {
103   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)(A)->data;
104   Mat_SuperLU    *lu = (Mat_SuperLU*)(*F)->spptr;
105   PetscErrorCode ierr;
106   PetscInt       sinfo;
107   SuperLUStat_t  stat;
108   PetscReal      ferr, berr;
109   NCformat       *Ustore;
110   SCformat       *Lstore;
111 
112   PetscFunctionBegin;
113   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
114     lu->options.Fact = SamePattern;
115     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
116     Destroy_SuperMatrix_Store(&lu->A);
117     if ( lu->lwork >= 0 ) {
118       Destroy_SuperNode_Matrix(&lu->L);
119       Destroy_CompCol_Matrix(&lu->U);
120       lu->options.Fact = SamePattern;
121     }
122   }
123 
124   /* Create the SuperMatrix for lu->A=A^T:
125        Since SuperLU likes column-oriented matrices,we pass it the transpose,
126        and then solve A^T X = B in MatSolve(). */
127 #if defined(PETSC_USE_COMPLEX)
128   zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
129                            SLU_NC,SLU_Z,SLU_GE);
130 #else
131   dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
132                            SLU_NC,SLU_D,SLU_GE);
133 #endif
134 
135   /* Initialize the statistics variables. */
136   StatInit(&stat);
137 
138   /* Numerical factorization */
139   lu->B.ncol = 0;  /* Indicate not to solve the system */
140 #if defined(PETSC_USE_COMPLEX)
141    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
142            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
143            &lu->mem_usage, &stat, &sinfo);
144 #else
145   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
146            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
147            &lu->mem_usage, &stat, &sinfo);
148 #endif
149   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
150     if ( lu->options.PivotGrowth )
151       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
152     if ( lu->options.ConditionNumber )
153       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
154   } else if ( sinfo > 0 ){
155     if ( lu->lwork == -1 ) {
156       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
157     } else {
158       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
159     }
160   } else { /* sinfo < 0 */
161     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
162   }
163 
164   if ( lu->options.PrintStat ) {
165     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
166     StatPrint(&stat);
167     Lstore = (SCformat *) lu->L.Store;
168     Ustore = (NCformat *) lu->U.Store;
169     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
170     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
171     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
172     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n",
173 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
174 	       lu->mem_usage.expansions);
175   }
176   StatFree(&stat);
177 
178   lu->flg = SAME_NONZERO_PATTERN;
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "MatDestroy_SuperLU"
184 PetscErrorCode MatDestroy_SuperLU(Mat A)
185 {
186   PetscErrorCode ierr;
187   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
188 
189   PetscFunctionBegin;
190   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
191     Destroy_SuperMatrix_Store(&lu->A);
192     Destroy_SuperMatrix_Store(&lu->B);
193     Destroy_SuperMatrix_Store(&lu->X);
194 
195     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
196     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
197     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
198     ierr = PetscFree(lu->R);CHKERRQ(ierr);
199     ierr = PetscFree(lu->C);CHKERRQ(ierr);
200     if ( lu->lwork >= 0 ) {
201       Destroy_SuperNode_Matrix(&lu->L);
202       Destroy_CompCol_Matrix(&lu->U);
203     }
204   }
205   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "MatView_SuperLU"
211 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
212 {
213   PetscErrorCode    ierr;
214   PetscTruth        iascii;
215   PetscViewerFormat format;
216 
217   PetscFunctionBegin;
218   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
219 
220   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
221   if (iascii) {
222     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
223     if (format == PETSC_VIEWER_ASCII_INFO) {
224       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
225     }
226   }
227   PetscFunctionReturn(0);
228 }
229 
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "MatSolve_SuperLU_Private"
233 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
234 {
235   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
236   PetscScalar    *barray,*xarray;
237   PetscErrorCode ierr;
238   PetscInt       info,i;
239   SuperLUStat_t  stat;
240   PetscReal      ferr,berr;
241 
242   PetscFunctionBegin;
243   if ( lu->lwork == -1 ) {
244     PetscFunctionReturn(0);
245   }
246   lu->B.ncol = 1;   /* Set the number of right-hand side */
247   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
248   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
249 
250 #if defined(PETSC_USE_COMPLEX)
251   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
252   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
253 #else
254   ((DNformat*)lu->B.Store)->nzval = barray;
255   ((DNformat*)lu->X.Store)->nzval = xarray;
256 #endif
257 
258   /* Initialize the statistics variables. */
259   StatInit(&stat);
260 
261   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
262 #if defined(PETSC_USE_COMPLEX)
263   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
264            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
265            &lu->mem_usage, &stat, &info);
266 #else
267   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
268            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
269            &lu->mem_usage, &stat, &info);
270 #endif
271   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
272   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
273 
274   if ( !info || info == lu->A.ncol+1 ) {
275     if ( lu->options.IterRefine ) {
276       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
277       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
278       for (i = 0; i < 1; ++i)
279         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
280     }
281   } else if ( info > 0 ){
282     if ( lu->lwork == -1 ) {
283       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
284     } else {
285       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
286     }
287   } else if (info < 0){
288     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
289   }
290 
291   if ( lu->options.PrintStat ) {
292     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
293     StatPrint(&stat);
294   }
295   StatFree(&stat);
296   PetscFunctionReturn(0);
297 }
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "MatSolve_SuperLU"
301 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
302 {
303   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
304   PetscErrorCode ierr;
305 
306   PetscFunctionBegin;
307   lu->options.Trans = TRANS;
308   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "MatSolveTranspose_SuperLU"
314 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
315 {
316   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   lu->options.Trans = NOTRANS;
321   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
322   PetscFunctionReturn(0);
323 }
324 
325 
326 /*
327    Note the r permutation is ignored
328 */
329 #undef __FUNCT__
330 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
331 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
332 {
333   Mat_SuperLU    *lu = (Mat_SuperLU*)((*F)->spptr);
334   PetscErrorCode ierr;
335   PetscInt       m=A->rmap->n,n=A->cmap->n;
336 
337   PetscFunctionBegin;
338 
339   /* Allocate spaces (notice sizes are for the transpose) */
340   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
341   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
342   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
343   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr);
344   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr);
345 
346   /* create rhs and solution x without allocate space for .Store */
347 #if defined(PETSC_USE_COMPLEX)
348   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
349   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
350 #else
351   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
352   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
353 #endif
354 
355   lu->flg            = DIFFERENT_NONZERO_PATTERN;
356   lu->CleanUpSuperLU = PETSC_TRUE;
357   (*F)->ops->lufactornumeric  = MatLUFactorNumeric_SuperLU;
358   (*F)->ops->solve            = MatSolve_SuperLU;
359   (*F)->ops->solvetranspose   = MatSolveTranspose_SuperLU;
360   PetscFunctionReturn(0);
361 }
362 
363 
364 /*MC
365   MAT_SOLVER_SUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
366   via the external package SuperLU.
367 
368   Use config/configure.py --download-superlu to have PETSc installed with SuperLU
369 
370   Options Database Keys:
371 + -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
372                                     1: MMD applied to A'*A,
373                                     2: MMD applied to A'+A,
374                                     3: COLAMD, approximate minimum degree column ordering
375 . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
376                           choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
377 - -mat_superlu_printstat - print SuperLU statistics about the factorization
378 
379    Notes: Do not confuse this with MAT_SOLVER_SUPERLU_DIST which is for parallel sparse solves
380 
381    Level: beginner
382 
383 .seealso: PCLU, MAT_SOLVER_SUPERLU_DIST, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES
384 M*/
385 
386 EXTERN_C_BEGIN
387 #undef __FUNCT__
388 #define __FUNCT__ "MatGetFactor_seqaij_superlu"
389 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
390 {
391   Mat            B;
392   Mat_SuperLU    *lu;
393   PetscErrorCode ierr;
394   PetscInt       indx;
395   PetscTruth     flg;
396   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
397   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
398   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
399 
400   PetscFunctionBegin;
401   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
402   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
403   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
404   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
405 
406   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
407   B->ops->destroy          = MatDestroy_SuperLU;
408   B->factor               = MAT_FACTOR_LU;
409   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
410   B->preallocated         = PETSC_TRUE;
411 
412   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
413   set_default_options(&lu->options);
414   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
415   lu->options.Equil = NO;
416   lu->options.PrintStat = NO;
417   lu->lwork = 0;   /* allocate space internally by system malloc */
418 
419   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
420     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
421     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
422     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
423     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
424     ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
425     if (flg) lu->options.SymmetricMode = YES;
426     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
427     ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
428     if (flg) lu->options.PivotGrowth = YES;
429     ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
430     if (flg) lu->options.ConditionNumber = YES;
431     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
432     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
433     ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
434     if (flg) lu->options.ReplaceTinyPivot = YES;
435     ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
436     if (flg) lu->options.PrintStat = YES;
437     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
438     if (lu->lwork > 0 ){
439       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
440     } else if (lu->lwork != 0 && lu->lwork != -1){
441       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
442       lu->lwork = 0;
443     }
444   PetscOptionsEnd();
445 
446 #ifdef SUPERLU2
447   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
448 #endif
449   B->spptr = lu;
450   *F = B;
451   PetscFunctionReturn(0);
452 }
453 EXTERN_C_END
454