xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision f3fe499b4cc4d64bf04aa4f5e4963dcc4eb56541)
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   PetscTruth 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 
225     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
226     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
227     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
228     ierr = PetscFree(lu->R);CHKERRQ(ierr);
229     ierr = PetscFree(lu->C);CHKERRQ(ierr);
230     if ( lu->lwork >= 0 ) {
231       Destroy_SuperNode_Matrix(&lu->L);
232       Destroy_CompCol_Matrix(&lu->U);
233     }
234   }
235   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
236   if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);}
237   if (lu->rhs_dup){ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);}
238   PetscFunctionReturn(0);
239 }
240 
241 #undef __FUNCT__
242 #define __FUNCT__ "MatView_SuperLU"
243 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
244 {
245   PetscErrorCode    ierr;
246   PetscTruth        iascii;
247   PetscViewerFormat format;
248 
249   PetscFunctionBegin;
250   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
251 
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__ "MatFactorGetSolverPackage_seqaij_superlu"
400 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
401 {
402   PetscFunctionBegin;
403   *type = MATSOLVERSUPERLU;
404   PetscFunctionReturn(0);
405 }
406 EXTERN_C_END
407 
408 
409 /*MC
410   MATSOLVERSUPERLU = "superlu" - A solver package roviding direct solvers (LU) for sequential matrices
411   via the external package SuperLU.
412 
413   Use ./configure --download-superlu to have PETSc installed with SuperLU
414 
415   Options Database Keys:
416 + -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
417                                     1: MMD applied to A'*A,
418                                     2: MMD applied to A'+A,
419                                     3: COLAMD, approximate minimum degree column ordering
420 . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
421                           choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
422 - -mat_superlu_printstat - print SuperLU statistics about the factorization
423 
424    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
425 
426    Level: beginner
427 
428 .seealso: PCLU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
429 M*/
430 
431 EXTERN_C_BEGIN
432 #undef __FUNCT__
433 #define __FUNCT__ "MatGetFactor_seqaij_superlu"
434 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
435 {
436   Mat            B;
437   Mat_SuperLU    *lu;
438   PetscErrorCode ierr;
439   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
440   PetscTruth     flg;
441   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
442   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
443   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
444 
445   PetscFunctionBegin;
446   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
447   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
448   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
449   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
450 
451   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
452     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
453     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
454   } else {
455     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
456   }
457 
458   B->ops->destroy          = MatDestroy_SuperLU;
459   B->ops->view             = MatView_SuperLU;
460   B->factortype            = ftype;
461   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
462   B->preallocated          = PETSC_TRUE;
463 
464   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
465 
466   if (ftype == MAT_FACTOR_LU){
467     set_default_options(&lu->options);
468     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
469       "Whether or not the system will be equilibrated depends on the
470        scaling of the matrix A, but if equilibration is used, A is
471        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
472        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
473      We set 'options.Equil = NO' as default because additional space is needed for it.
474     */
475     lu->options.Equil     = NO;
476   } else if (ftype == MAT_FACTOR_ILU){
477     /* Set the default input options of ilu:
478 	options.Fact = DOFACT;
479 	options.Equil = YES;           // must be YES for ilu - don't know why
480 	options.ColPerm = COLAMD;
481 	options.DiagPivotThresh = 0.1; //different from complete LU
482 	options.Trans = NOTRANS;
483 	options.IterRefine = NOREFINE;
484 	options.SymmetricMode = NO;
485 	options.PivotGrowth = NO;
486 	options.ConditionNumber = NO;
487 	options.PrintStat = YES;
488 	options.RowPerm = LargeDiag;
489 	options.ILU_DropTol = 1e-4;
490 	options.ILU_FillTol = 1e-2;
491 	options.ILU_FillFactor = 10.0;
492 	options.ILU_DropRule = DROP_BASIC | DROP_AREA;
493 	options.ILU_Norm = INF_NORM;
494 	options.ILU_MILU = SMILU_2;
495     */
496     ilu_set_default_options(&lu->options);
497   }
498   lu->options.PrintStat = NO;
499 
500   /* Initialize the statistics variables. */
501   StatInit(&lu->stat);
502   lu->lwork = 0;   /* allocate space internally by system malloc */
503 
504   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
505     ierr = PetscOptionsTruth("-mat_superlu_equil","Equil","None",(PetscTruth)lu->options.Equil,(PetscTruth*)&lu->options.Equil,0);CHKERRQ(ierr);
506     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
507     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
508     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
509     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
510     ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscTruth)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
511     if (flg) lu->options.SymmetricMode = YES;
512     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
513     ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscTruth)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
514     if (flg) lu->options.PivotGrowth = YES;
515     ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscTruth)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
516     if (flg) lu->options.ConditionNumber = YES;
517     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
518     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
519     ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscTruth)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
520     if (flg) lu->options.ReplaceTinyPivot = YES;
521     ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",(PetscTruth)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
522     if (flg) lu->options.PrintStat = YES;
523     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
524     if (lu->lwork > 0 ){
525       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
526     } else if (lu->lwork != 0 && lu->lwork != -1){
527       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
528       lu->lwork = 0;
529     }
530     /* ilu options */
531     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr);
532     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr);
533     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr);
534     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
535     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
536     if (flg){
537       lu->options.ILU_Norm = (norm_t)indx;
538     }
539     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
540     if (flg){
541       lu->options.ILU_MILU = (milu_t)indx;
542     }
543   PetscOptionsEnd();
544   if (lu->options.Equil == YES) {
545     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
546     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
547   }
548 
549   /* Allocate spaces (notice sizes are for the transpose) */
550   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
551   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
552   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
553   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
554   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
555 
556   /* create rhs and solution x without allocate space for .Store */
557 #if defined(PETSC_USE_COMPLEX)
558   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
559   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
560 #else
561   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
562   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
563 #endif
564 
565 #ifdef SUPERLU2
566   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
567 #endif
568   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
569   B->spptr = lu;
570   *F = B;
571   PetscFunctionReturn(0);
572 }
573 EXTERN_C_END
574