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