xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 3c48a1e8da19189ff2402a4e41a2fc082d52c349)
1 
2 /*  --------------------------------------------------------------------
3 
4      This file implements a subclass of the SeqAIJ matrix class that uses
5      the SuperLU sparse solver. You can use this as a starting point for
6      implementing your own subclass of a PETSc matrix class.
7 
8      This demonstrates a way to make an implementation inheritence of a PETSc
9      matrix type. This means constructing a new matrix type (SuperLU) by changing some
10      of the methods of the previous type (SeqAIJ), adding additional data, and possibly
11      additional method. (See any book on object oriented programming).
12 */
13 
14 /*
15      Defines the data structure for the base matrix type (SeqAIJ)
16 */
17 #include "../src/mat/impls/aij/seq/aij.h"    /*I "petscmat.h" I*/
18 
19 /*
20      SuperLU include files
21 */
22 EXTERN_C_BEGIN
23 #if defined(PETSC_USE_COMPLEX)
24 #include "slu_zdefs.h"
25 #else
26 #include "slu_ddefs.h"
27 #endif
28 #include "slu_util.h"
29 EXTERN_C_END
30 
31 /*
32      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
33 */
34 typedef struct {
35   SuperMatrix       A,L,U,B,X;
36   superlu_options_t options;
37   PetscInt          *perm_c; /* column permutation vector */
38   PetscInt          *perm_r; /* row permutations from partial pivoting */
39   PetscInt          *etree;
40   PetscReal         *R, *C;
41   char              equed[1];
42   PetscInt          lwork;
43   void              *work;
44   PetscReal         rpg, rcond;
45   mem_usage_t       mem_usage;
46   MatStructure      flg;
47   SuperLUStat_t     stat;
48   Mat               A_dup;
49   PetscScalar       *rhs_dup;
50 
51   /* Flag to clean up (non-global) SuperLU objects during Destroy */
52   PetscBool  CleanUpSuperLU;
53 } Mat_SuperLU;
54 
55 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
56 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *);
57 extern PetscErrorCode MatDestroy_SuperLU(Mat);
58 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
59 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
60 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
61 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
62 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*);
63 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *);
64 
65 /*
66     Utility function
67 */
68 #undef __FUNCT__
69 #define __FUNCT__ "MatFactorInfo_SuperLU"
70 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
71 {
72   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
73   PetscErrorCode    ierr;
74   superlu_options_t options;
75 
76   PetscFunctionBegin;
77   /* check if matrix is superlu_dist type */
78   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
79 
80   options = lu->options;
81   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
82   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
83   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
84   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
85   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
86   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
87   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
88   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
89   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
90   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
91   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
92   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
93   if (A->factortype == MAT_FACTOR_ILU){
94     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr);
95     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr);
96     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr);
97     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr);
98     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr);
99     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr);
100   }
101   PetscFunctionReturn(0);
102 }
103 
104 /*
105     These are the methods provided to REPLACE the corresponding methods of the
106    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
107 */
108 #undef __FUNCT__
109 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
110 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
111 {
112   Mat_SuperLU    *lu = (Mat_SuperLU*)F->spptr;
113   Mat_SeqAIJ     *aa;
114   PetscErrorCode ierr;
115   PetscInt       sinfo;
116   PetscReal      ferr, berr;
117   NCformat       *Ustore;
118   SCformat       *Lstore;
119 
120   PetscFunctionBegin;
121   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
122     lu->options.Fact = SamePattern;
123     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
124     Destroy_SuperMatrix_Store(&lu->A);
125     if (lu->options.Equil){
126       ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
127     }
128     if ( lu->lwork >= 0 ) {
129       Destroy_SuperNode_Matrix(&lu->L);
130       Destroy_CompCol_Matrix(&lu->U);
131       lu->options.Fact = SamePattern;
132     }
133   }
134 
135   /* Create the SuperMatrix for lu->A=A^T:
136        Since SuperLU likes column-oriented matrices,we pass it the transpose,
137        and then solve A^T X = B in MatSolve(). */
138   if (lu->options.Equil){
139     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
140   } else {
141     aa = (Mat_SeqAIJ*)(A)->data;
142   }
143 #if defined(PETSC_USE_COMPLEX)
144   zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
145                            SLU_NC,SLU_Z,SLU_GE);
146 #else
147   dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
148                            SLU_NC,SLU_D,SLU_GE);
149 #endif
150 
151   /* Numerical factorization */
152   lu->B.ncol = 0;  /* Indicate not to solve the system */
153   if (F->factortype == MAT_FACTOR_LU){
154 #if defined(PETSC_USE_COMPLEX)
155     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
156            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
157            &lu->mem_usage, &lu->stat, &sinfo);
158 #else
159     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
160            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
161            &lu->mem_usage, &lu->stat, &sinfo);
162 #endif
163   } else if (F->factortype == MAT_FACTOR_ILU){
164     /* Compute the incomplete factorization, condition number and pivot growth */
165 #if defined(PETSC_USE_COMPLEX)
166     zgsisx(&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,
168            &lu->mem_usage, &lu->stat, &sinfo);
169 #else
170     dgsisx(&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,
172           &lu->mem_usage, &lu->stat, &sinfo);
173 #endif
174   } else {
175     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
176   }
177   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
178     if ( lu->options.PivotGrowth )
179       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
180     if ( lu->options.ConditionNumber )
181       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
182   } else if ( sinfo > 0 ){
183     if ( lu->lwork == -1 ) {
184       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
185     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
186   } else { /* sinfo < 0 */
187     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
188   }
189 
190   if ( lu->options.PrintStat ) {
191     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
192     StatPrint(&lu->stat);
193     Lstore = (SCformat *) lu->L.Store;
194     Ustore = (NCformat *) lu->U.Store;
195     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
196     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
197     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
198     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
199 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
200   }
201 
202   lu->flg = SAME_NONZERO_PATTERN;
203   F->ops->solve          = MatSolve_SuperLU;
204   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "MatDestroy_SuperLU"
210 PetscErrorCode MatDestroy_SuperLU(Mat A)
211 {
212   PetscErrorCode ierr;
213   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
214 
215   PetscFunctionBegin;
216   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
217     Destroy_SuperMatrix_Store(&lu->A);
218     Destroy_SuperMatrix_Store(&lu->B);
219     Destroy_SuperMatrix_Store(&lu->X);
220     StatFree(&lu->stat);
221     if ( lu->lwork >= 0 ) {
222       Destroy_SuperNode_Matrix(&lu->L);
223       Destroy_CompCol_Matrix(&lu->U);
224     }
225   }
226 
227   ierr = PetscFree(lu->etree);CHKERRQ(ierr);
228   ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
229   ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
230   ierr = PetscFree(lu->R);CHKERRQ(ierr);
231   ierr = PetscFree(lu->C);CHKERRQ(ierr);
232 
233   /* clear composed functions */
234   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
235   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr);
236 
237   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
238   if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);}
239   ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "MatView_SuperLU"
245 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
246 {
247   PetscErrorCode    ierr;
248   PetscBool         iascii;
249   PetscViewerFormat format;
250 
251   PetscFunctionBegin;
252   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
253   if (iascii) {
254     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
255     if (format == PETSC_VIEWER_ASCII_INFO) {
256       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
257     }
258   }
259   PetscFunctionReturn(0);
260 }
261 
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "MatSolve_SuperLU_Private"
265 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
266 {
267   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
268   PetscScalar    *barray,*xarray;
269   PetscErrorCode ierr;
270   PetscInt       info,i,n=x->map->n;
271   PetscReal      ferr,berr;
272 
273   PetscFunctionBegin;
274   if ( lu->lwork == -1 ) {
275     PetscFunctionReturn(0);
276   }
277 
278   lu->B.ncol = 1;   /* Set the number of right-hand side */
279   if (lu->options.Equil && !lu->rhs_dup){
280     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
281     ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr);
282   }
283   if (lu->options.Equil){
284     /* Copy b into rsh_dup */
285     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
286     ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
287     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
288     barray = lu->rhs_dup;
289   } else {
290     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
291   }
292   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
293 
294 #if defined(PETSC_USE_COMPLEX)
295   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
296   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
297 #else
298   ((DNformat*)lu->B.Store)->nzval = barray;
299   ((DNformat*)lu->X.Store)->nzval = xarray;
300 #endif
301 
302   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
303   if (A->factortype == MAT_FACTOR_LU){
304 #if defined(PETSC_USE_COMPLEX)
305     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
306            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
307            &lu->mem_usage, &lu->stat, &info);
308 #else
309     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
310            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
311            &lu->mem_usage, &lu->stat, &info);
312 #endif
313   } else if (A->factortype == MAT_FACTOR_ILU){
314 #if defined(PETSC_USE_COMPLEX)
315     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
316            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
317            &lu->mem_usage, &lu->stat, &info);
318 #else
319     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
320            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
321            &lu->mem_usage, &lu->stat, &info);
322 #endif
323   } else {
324     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
325   }
326   if (!lu->options.Equil){
327     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
328   }
329   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
330 
331   if ( !info || info == lu->A.ncol+1 ) {
332     if ( lu->options.IterRefine ) {
333       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
334       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
335       for (i = 0; i < 1; ++i)
336         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
337     }
338   } else if ( info > 0 ){
339     if ( lu->lwork == -1 ) {
340       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
341     } else {
342       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
343     }
344   } else if (info < 0){
345     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
346   }
347 
348   if ( lu->options.PrintStat ) {
349     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
350     StatPrint(&lu->stat);
351   }
352   PetscFunctionReturn(0);
353 }
354 
355 #undef __FUNCT__
356 #define __FUNCT__ "MatSolve_SuperLU"
357 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
358 {
359   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
360   PetscErrorCode ierr;
361 
362   PetscFunctionBegin;
363   lu->options.Trans = TRANS;
364   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
365   PetscFunctionReturn(0);
366 }
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "MatSolveTranspose_SuperLU"
370 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
371 {
372   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
373   PetscErrorCode ierr;
374 
375   PetscFunctionBegin;
376   lu->options.Trans = NOTRANS;
377   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 /*
382    Note the r permutation is ignored
383 */
384 #undef __FUNCT__
385 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
386 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
387 {
388   Mat_SuperLU    *lu = (Mat_SuperLU*)(F->spptr);
389 
390   PetscFunctionBegin;
391   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
392   lu->CleanUpSuperLU      = PETSC_TRUE;
393   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
394   PetscFunctionReturn(0);
395 }
396 
397 EXTERN_C_BEGIN
398 #undef __FUNCT__
399 #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
400 PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
401 {
402   Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
403 
404   PetscFunctionBegin;
405   lu->options.ILU_DropTol = dtol;
406   PetscFunctionReturn(0);
407 }
408 EXTERN_C_END
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "MatSuperluSetILUDropTol"
412 /*@
413   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
414    Logically Collective on Mat
415 
416    Input Parameters:
417 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
418 -  dtol - drop tolerance
419 
420   Options Database:
421 .   -mat_superlu_ilu_droptol <dtol>
422 
423    Level: beginner
424 
425    References: SuperLU Users' Guide
426 
427 .seealso: MatGetFactor()
428 @*/
429 PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
430 {
431   PetscErrorCode ierr;
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
435   PetscValidLogicalCollectiveInt(F,dtol,2);
436   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
437   PetscFunctionReturn(0);
438 }
439 
440 EXTERN_C_BEGIN
441 #undef __FUNCT__
442 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
443 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
444 {
445   PetscFunctionBegin;
446   *type = MATSOLVERSUPERLU;
447   PetscFunctionReturn(0);
448 }
449 EXTERN_C_END
450 
451 
452 /*MC
453   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
454   via the external package SuperLU.
455 
456   Use ./configure --download-superlu to have PETSc installed with SuperLU
457 
458   Options Database Keys:
459 +  -mat_superlu_equil: <FALSE> Equil (None)
460 .  -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
461 .  -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA
462 .  -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None)
463 .  -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None)
464 .  -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None)
465 .  -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None)
466 .  -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag
467 .  -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None)
468 .  -mat_superlu_printstat: <FALSE> PrintStat (None)
469 .  -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None)
470 .  -mat_superlu_ilu_droptol <0>: ILU_DropTol (None)
471 .  -mat_superlu_ilu_filltol <0>: ILU_FillTol (None)
472 .  -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None)
473 .  -mat_superlu_ilu_droprull <0>: ILU_DropRule (None)
474 .  -mat_superlu_ilu_norm <0>: ILU_Norm (None)
475 -  -mat_superlu_ilu_milu <0>: ILU_MILU (None)
476 
477    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
478 
479    Level: beginner
480 
481 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
482 M*/
483 
484 EXTERN_C_BEGIN
485 #undef __FUNCT__
486 #define __FUNCT__ "MatGetFactor_seqaij_superlu"
487 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
488 {
489   Mat            B;
490   Mat_SuperLU    *lu;
491   PetscErrorCode ierr;
492   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
493   PetscBool      flg;
494   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
495   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
496   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
497 
498   PetscFunctionBegin;
499   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
500   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
501   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
502   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
503 
504   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
505     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
506     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
507   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
508 
509   B->ops->destroy          = MatDestroy_SuperLU;
510   B->ops->view             = MatView_SuperLU;
511   B->factortype            = ftype;
512   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
513   B->preallocated          = PETSC_TRUE;
514 
515   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
516 
517   if (ftype == MAT_FACTOR_LU){
518     set_default_options(&lu->options);
519     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
520       "Whether or not the system will be equilibrated depends on the
521        scaling of the matrix A, but if equilibration is used, A is
522        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
523        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
524      We set 'options.Equil = NO' as default because additional space is needed for it.
525     */
526     lu->options.Equil = NO;
527   } else if (ftype == MAT_FACTOR_ILU){
528     /* Set the default input options of ilu:
529 	options.Fact = DOFACT;
530 	options.Equil = YES;           // must be YES for ilu - don't know why
531 	options.ColPerm = COLAMD;
532 	options.DiagPivotThresh = 0.1; //different from complete LU
533 	options.Trans = NOTRANS;
534 	options.IterRefine = NOREFINE;
535 	options.SymmetricMode = NO;
536 	options.PivotGrowth = NO;
537 	options.ConditionNumber = NO;
538 	options.PrintStat = YES;
539 	options.RowPerm = LargeDiag;
540 	options.ILU_DropTol = 1e-4;
541 	options.ILU_FillTol = 1e-2;
542 	options.ILU_FillFactor = 10.0;
543 	options.ILU_DropRule = DROP_BASIC | DROP_AREA;
544 	options.ILU_Norm = INF_NORM;
545 	options.ILU_MILU = SMILU_2;
546     */
547     ilu_set_default_options(&lu->options);
548   }
549   lu->options.PrintStat = NO;
550 
551   /* Initialize the statistics variables. */
552   StatInit(&lu->stat);
553   lu->lwork = 0;   /* allocate space internally by system malloc */
554 
555   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
556     ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr);
557     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
558     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
559     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
560     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
561     ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
562     if (flg) lu->options.SymmetricMode = YES;
563     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
564     ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
565     if (flg) lu->options.PivotGrowth = YES;
566     ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
567     if (flg) lu->options.ConditionNumber = YES;
568     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
569     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
570     ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
571     if (flg) lu->options.ReplaceTinyPivot = YES;
572     ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
573     if (flg) lu->options.PrintStat = YES;
574     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
575     if (lu->lwork > 0 ){
576       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
577     } else if (lu->lwork != 0 && lu->lwork != -1){
578       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
579       lu->lwork = 0;
580     }
581     /* ilu options */
582     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr);
583     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr);
584     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr);
585     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
586     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
587     if (flg){
588       lu->options.ILU_Norm = (norm_t)indx;
589     }
590     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
591     if (flg){
592       lu->options.ILU_MILU = (milu_t)indx;
593     }
594   PetscOptionsEnd();
595   if (lu->options.Equil == YES) {
596     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
597     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
598   }
599 
600   /* Allocate spaces (notice sizes are for the transpose) */
601   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
602   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
603   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
604   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
605   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
606 
607   /* create rhs and solution x without allocate space for .Store */
608 #if defined(PETSC_USE_COMPLEX)
609   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
610   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
611 #else
612   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
613   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
614 #endif
615 
616 #ifdef SUPERLU2
617   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
618 #endif
619   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
620   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
621   B->spptr = lu;
622   *F = B;
623   PetscFunctionReturn(0);
624 }
625 EXTERN_C_END
626 
627