xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision efe48dd8cf8b935accbbb9f4bcb20bc83865fa4d)
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"
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 {
187       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
188     }
189   } else { /* sinfo < 0 */
190     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
191   }
192 
193   if ( lu->options.PrintStat ) {
194     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
195     StatPrint(&lu->stat);
196     Lstore = (SCformat *) lu->L.Store;
197     Ustore = (NCformat *) lu->U.Store;
198     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
199     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
200     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
201     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
202 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
203   }
204 
205   lu->flg = SAME_NONZERO_PATTERN;
206   (F)->ops->solve          = MatSolve_SuperLU;
207   (F)->ops->solvetranspose = MatSolveTranspose_SuperLU;
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatDestroy_SuperLU"
213 PetscErrorCode MatDestroy_SuperLU(Mat A)
214 {
215   PetscErrorCode ierr;
216   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
217 
218   PetscFunctionBegin;
219   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
220     Destroy_SuperMatrix_Store(&lu->A);
221     Destroy_SuperMatrix_Store(&lu->B);
222     Destroy_SuperMatrix_Store(&lu->X);
223     StatFree(&lu->stat);
224     if ( lu->lwork >= 0 ) {
225       Destroy_SuperNode_Matrix(&lu->L);
226       Destroy_CompCol_Matrix(&lu->U);
227     }
228 
229   }
230 
231   ierr = PetscFree(lu->etree);CHKERRQ(ierr);
232   ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
233   ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
234   ierr = PetscFree(lu->R);CHKERRQ(ierr);
235   ierr = PetscFree(lu->C);CHKERRQ(ierr);
236 
237   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
238   if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);}
239   if (lu->rhs_dup){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 = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
253 
254   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
255   if (iascii) {
256     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
257     if (format == PETSC_VIEWER_ASCII_INFO) {
258       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
259     }
260   }
261   PetscFunctionReturn(0);
262 }
263 
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "MatSolve_SuperLU_Private"
267 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
268 {
269   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
270   PetscScalar    *barray,*xarray;
271   PetscErrorCode ierr;
272   PetscInt       info,i,n=x->map->n;
273   PetscReal      ferr,berr;
274 
275   PetscFunctionBegin;
276   if ( lu->lwork == -1 ) {
277     PetscFunctionReturn(0);
278   }
279 
280   lu->B.ncol = 1;   /* Set the number of right-hand side */
281   if (lu->options.Equil && !lu->rhs_dup){
282     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
283     ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr);
284   }
285   if (lu->options.Equil){
286     /* Copy b into rsh_dup */
287     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
288     ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
289     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
290     barray = lu->rhs_dup;
291   } else {
292     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
293   }
294   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
295 
296 #if defined(PETSC_USE_COMPLEX)
297   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
298   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
299 #else
300   ((DNformat*)lu->B.Store)->nzval = barray;
301   ((DNformat*)lu->X.Store)->nzval = xarray;
302 #endif
303 
304   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
305   if (A->factortype == MAT_FACTOR_LU){
306 #if defined(PETSC_USE_COMPLEX)
307     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
308            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
309            &lu->mem_usage, &lu->stat, &info);
310 #else
311     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
312            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
313            &lu->mem_usage, &lu->stat, &info);
314 #endif
315   } else if (A->factortype == MAT_FACTOR_ILU){
316 #if defined(PETSC_USE_COMPLEX)
317     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
318            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
319            &lu->mem_usage, &lu->stat, &info);
320 #else
321     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
322            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
323            &lu->mem_usage, &lu->stat, &info);
324 #endif
325   } else {
326     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
327   }
328   if (!lu->options.Equil){
329     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
330   }
331   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
332 
333   if ( !info || info == lu->A.ncol+1 ) {
334     if ( lu->options.IterRefine ) {
335       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
336       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
337       for (i = 0; i < 1; ++i)
338         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
339     }
340   } else if ( info > 0 ){
341     if ( lu->lwork == -1 ) {
342       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
343     } else {
344       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
345     }
346   } else if (info < 0){
347     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
348   }
349 
350   if ( lu->options.PrintStat ) {
351     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
352     StatPrint(&lu->stat);
353   }
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "MatSolve_SuperLU"
359 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
360 {
361   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   lu->options.Trans = TRANS;
366   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "MatSolveTranspose_SuperLU"
372 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
373 {
374   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
375   PetscErrorCode ierr;
376 
377   PetscFunctionBegin;
378   lu->options.Trans = NOTRANS;
379   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 /*
384    Note the r permutation is ignored
385 */
386 #undef __FUNCT__
387 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
388 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
389 {
390   Mat_SuperLU    *lu = (Mat_SuperLU*)((F)->spptr);
391 
392   PetscFunctionBegin;
393   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
394   lu->CleanUpSuperLU        = PETSC_TRUE;
395   (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
396   PetscFunctionReturn(0);
397 }
398 
399 EXTERN_C_BEGIN
400 #undef __FUNCT__
401 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
402 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
403 {
404   PetscFunctionBegin;
405   *type = MATSOLVERSUPERLU;
406   PetscFunctionReturn(0);
407 }
408 EXTERN_C_END
409 
410 
411 /*MC
412   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
413   via the external package SuperLU.
414 
415   Use ./configure --download-superlu to have PETSc installed with SuperLU
416 
417   Options Database Keys:
418 +  -mat_superlu_equil: <FALSE> Equil (None)
419 .  -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
420 .  -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA
421 .  -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None)
422 .  -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None)
423 .  -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None)
424 .  -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None)
425 .  -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag
426 .  -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None)
427 .  -mat_superlu_printstat: <FALSE> PrintStat (None)
428 .  -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None)
429 .  -mat_superlu_ilu_droptol <0>: ILU_DropTol (None)
430 .  -mat_superlu_ilu_filltol <0>: ILU_FillTol (None)
431 .  -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None)
432 .  -mat_superlu_ilu_droprull <0>: ILU_DropRule (None)
433 .  -mat_superlu_ilu_norm <0>: ILU_Norm (None)
434 -  -mat_superlu_ilu_milu <0>: ILU_MILU (None)
435 
436    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
437 
438    Level: beginner
439 
440 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
441 M*/
442 
443 EXTERN_C_BEGIN
444 #undef __FUNCT__
445 #define __FUNCT__ "MatGetFactor_seqaij_superlu"
446 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
447 {
448   Mat            B;
449   Mat_SuperLU    *lu;
450   PetscErrorCode ierr;
451   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
452   PetscBool      flg;
453   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
454   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
455   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
456 
457   PetscFunctionBegin;
458   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
459   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
460   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
461   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
462 
463   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
464     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
465     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
466   } else {
467     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
468   }
469 
470   B->ops->destroy          = MatDestroy_SuperLU;
471   B->ops->view             = MatView_SuperLU;
472   B->factortype            = ftype;
473   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
474   B->preallocated          = PETSC_TRUE;
475 
476   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
477 
478   if (ftype == MAT_FACTOR_LU){
479     set_default_options(&lu->options);
480     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
481       "Whether or not the system will be equilibrated depends on the
482        scaling of the matrix A, but if equilibration is used, A is
483        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
484        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
485      We set 'options.Equil = NO' as default because additional space is needed for it.
486     */
487     lu->options.Equil     = NO;
488   } else if (ftype == MAT_FACTOR_ILU){
489     /* Set the default input options of ilu:
490 	options.Fact = DOFACT;
491 	options.Equil = YES;           // must be YES for ilu - don't know why
492 	options.ColPerm = COLAMD;
493 	options.DiagPivotThresh = 0.1; //different from complete LU
494 	options.Trans = NOTRANS;
495 	options.IterRefine = NOREFINE;
496 	options.SymmetricMode = NO;
497 	options.PivotGrowth = NO;
498 	options.ConditionNumber = NO;
499 	options.PrintStat = YES;
500 	options.RowPerm = LargeDiag;
501 	options.ILU_DropTol = 1e-4;
502 	options.ILU_FillTol = 1e-2;
503 	options.ILU_FillFactor = 10.0;
504 	options.ILU_DropRule = DROP_BASIC | DROP_AREA;
505 	options.ILU_Norm = INF_NORM;
506 	options.ILU_MILU = SMILU_2;
507     */
508     ilu_set_default_options(&lu->options);
509   }
510   lu->options.PrintStat = NO;
511 
512   /* Initialize the statistics variables. */
513   StatInit(&lu->stat);
514   lu->lwork = 0;   /* allocate space internally by system malloc */
515 
516   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
517     ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr);
518     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
519     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
520     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
521     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
522     ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
523     if (flg) lu->options.SymmetricMode = YES;
524     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
525     ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
526     if (flg) lu->options.PivotGrowth = YES;
527     ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
528     if (flg) lu->options.ConditionNumber = YES;
529     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
530     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
531     ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
532     if (flg) lu->options.ReplaceTinyPivot = YES;
533     ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
534     if (flg) lu->options.PrintStat = YES;
535     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
536     if (lu->lwork > 0 ){
537       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
538     } else if (lu->lwork != 0 && lu->lwork != -1){
539       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
540       lu->lwork = 0;
541     }
542     /* ilu options */
543     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr);
544     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr);
545     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr);
546     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
547     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
548     if (flg){
549       lu->options.ILU_Norm = (norm_t)indx;
550     }
551     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
552     if (flg){
553       lu->options.ILU_MILU = (milu_t)indx;
554     }
555   PetscOptionsEnd();
556   if (lu->options.Equil == YES) {
557     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
558     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
559   }
560 
561   /* Allocate spaces (notice sizes are for the transpose) */
562   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
563   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
564   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
565   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
566   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
567 
568   /* create rhs and solution x without allocate space for .Store */
569 #if defined(PETSC_USE_COMPLEX)
570   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
571   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
572 #else
573   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
574   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
575 #endif
576 
577 #ifdef SUPERLU2
578   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
579 #endif
580   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
581   B->spptr = lu;
582   *F = B;
583   PetscFunctionReturn(0);
584 }
585 EXTERN_C_END
586