xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision e5bd52469b00a8f2942d86b46a52285e84b32048)
1 #define PETSCMAT_DLL
2 
3 /*  --------------------------------------------------------------------
4 
5      This file implements a subclass of the SeqAIJ matrix class that uses
6      the SuperLU 3.0 sparse solver. You can use this as a starting point for
7      implementing your own subclass of a PETSc matrix class.
8 
9      This demonstrates a way to make an implementation inheritence of a PETSc
10      matrix type. This means constructing a new matrix type (SuperLU) by changing some
11      of the methods of the previous type (SeqAIJ), adding additional data, and possibly
12      additional method. (See any book on object oriented programming).
13 */
14 
15 /*
16      Defines the data structure for the base matrix type (SeqAIJ)
17 */
18 #include "src/mat/impls/aij/seq/aij.h"
19 
20 /*
21      SuperLU include files
22 */
23 EXTERN_C_BEGIN
24 #if defined(PETSC_USE_COMPLEX)
25 #include "slu_zdefs.h"
26 #else
27 #include "slu_ddefs.h"
28 #endif
29 #include "slu_util.h"
30 EXTERN_C_END
31 
32 /*
33      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
34 */
35 typedef struct {
36   SuperMatrix       A,L,U,B,X;
37   superlu_options_t options;
38   PetscInt          *perm_c; /* column permutation vector */
39   PetscInt          *perm_r; /* row permutations from partial pivoting */
40   PetscInt          *etree;
41   PetscReal         *R, *C;
42   char              equed[1];
43   PetscInt          lwork;
44   void              *work;
45   PetscReal         rpg, rcond;
46   mem_usage_t       mem_usage;
47   MatStructure      flg;
48 
49   /*
50      This is where the methods for the superclass (SeqAIJ) are kept while we
51      reset the pointers in the function table to the new (SuperLU) versions
52   */
53   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
54   PetscErrorCode (*MatView)(Mat,PetscViewer);
55   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
56   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
57   PetscErrorCode (*MatDestroy)(Mat);
58 
59   /* Flag to clean up (non-global) SuperLU objects during Destroy */
60   PetscTruth CleanUpSuperLU;
61 } Mat_SuperLU;
62 
63 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
64 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,MatFactorInfo *,Mat *);
65 extern PetscErrorCode MatDestroy_SuperLU(Mat);
66 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
67 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
68 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
69 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
70 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo *,Mat *);
71 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *);
72 
73 /*
74     Takes a SuperLU matrix (that is a SeqAIJ matrix with the additional SuperLU data-structures
75    and methods) and converts it back to a regular SeqAIJ matrix.
76 */
77 EXTERN_C_BEGIN
78 #undef __FUNCT__
79 #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
80 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
81 {
82   PetscErrorCode ierr;
83   Mat            B=*newmat;
84   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
85 
86   PetscFunctionBegin;
87   if (reuse == MAT_INITIAL_MATRIX) {
88     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
89   }
90   /* Reset the original SeqAIJ function pointers */
91   B->ops->duplicate        = lu->MatDuplicate;
92   B->ops->view             = lu->MatView;
93   B->ops->assemblyend      = lu->MatAssemblyEnd;
94   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
95   B->ops->destroy          = lu->MatDestroy;
96 
97   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr);
98   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
99 
100   /* change the type name back to its original value */
101   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
102   *newmat = B;
103   PetscFunctionReturn(0);
104 }
105 EXTERN_C_END
106 
107 EXTERN_C_BEGIN
108 #undef __FUNCT__
109 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
110 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat)
111 {
112   PetscErrorCode ierr;
113   Mat            B=*newmat;
114   Mat_SuperLU    *lu;
115 
116   PetscFunctionBegin;
117   if (reuse == MAT_INITIAL_MATRIX) {
118     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
119   }
120 
121   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
122   /* save the original SeqAIJ methods that we are changing */
123   lu->MatDuplicate         = A->ops->duplicate;
124   lu->MatView              = A->ops->view;
125   lu->MatAssemblyEnd       = A->ops->assemblyend;
126   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
127   lu->MatDestroy           = A->ops->destroy;
128   lu->CleanUpSuperLU       = PETSC_FALSE;
129 
130   /* add to the matrix the location for all the SuperLU data is to be stored */
131   B->spptr                 = (void*)lu;
132 
133   /* set the methods in the function table to the SuperLU versions */
134   B->ops->duplicate        = MatDuplicate_SuperLU;
135   B->ops->view             = MatView_SuperLU;
136   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
137   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
138   B->ops->choleskyfactorsymbolic = 0;
139   B->ops->destroy          = MatDestroy_SuperLU;
140 
141   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
142                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
143   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
144                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
145   ierr = PetscInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
146   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
147   *newmat = B;
148   PetscFunctionReturn(0);
149 }
150 EXTERN_C_END
151 
152 /*
153     Utility function
154 */
155 #undef __FUNCT__
156 #define __FUNCT__ "MatFactorInfo_SuperLU"
157 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
158 {
159   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
160   PetscErrorCode    ierr;
161   superlu_options_t options;
162 
163   PetscFunctionBegin;
164   /* check if matrix is superlu_dist type */
165   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
166 
167   options = lu->options;
168   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
169   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
170   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
171   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
172   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
173   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
174   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
175   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
176   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
177   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
178   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
179   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
180 
181   PetscFunctionReturn(0);
182 }
183 
184 /*
185     These are the methods provided to REPLACE the corresponding methods of the
186    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
187 */
188 #undef __FUNCT__
189 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
190 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F)
191 {
192   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)(A)->data;
193   Mat_SuperLU    *lu = (Mat_SuperLU*)(*F)->spptr;
194   PetscErrorCode ierr;
195   PetscInt       sinfo;
196   SuperLUStat_t  stat;
197   PetscReal      ferr, berr;
198   NCformat       *Ustore;
199   SCformat       *Lstore;
200 
201   PetscFunctionBegin;
202   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
203     lu->options.Fact = SamePattern;
204     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
205     Destroy_SuperMatrix_Store(&lu->A);
206     if ( lu->lwork >= 0 ) {
207       Destroy_SuperNode_Matrix(&lu->L);
208       Destroy_CompCol_Matrix(&lu->U);
209       lu->options.Fact = SamePattern;
210     }
211   }
212 
213   /* Create the SuperMatrix for lu->A=A^T:
214        Since SuperLU likes column-oriented matrices,we pass it the transpose,
215        and then solve A^T X = B in MatSolve(). */
216 #if defined(PETSC_USE_COMPLEX)
217   zCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
218                            SLU_NC,SLU_Z,SLU_GE);
219 #else
220   dCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,aa->a,aa->j,aa->i,
221                            SLU_NC,SLU_D,SLU_GE);
222 #endif
223 
224   /* Initialize the statistics variables. */
225   StatInit(&stat);
226 
227   /* Numerical factorization */
228   lu->B.ncol = 0;  /* Indicate not to solve the system */
229 #if defined(PETSC_USE_COMPLEX)
230    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
231            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
232            &lu->mem_usage, &stat, &sinfo);
233 #else
234   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
235            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
236            &lu->mem_usage, &stat, &sinfo);
237 #endif
238   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
239     if ( lu->options.PivotGrowth )
240       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
241     if ( lu->options.ConditionNumber )
242       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
243   } else if ( sinfo > 0 ){
244     if ( lu->lwork == -1 ) {
245       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
246     } else {
247       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
248     }
249   } else { /* sinfo < 0 */
250     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
251   }
252 
253   if ( lu->options.PrintStat ) {
254     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
255     StatPrint(&stat);
256     Lstore = (SCformat *) lu->L.Store;
257     Ustore = (NCformat *) lu->U.Store;
258     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
259     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
260     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
261     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n",
262 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
263 	       lu->mem_usage.expansions);
264   }
265   StatFree(&stat);
266 
267   lu->flg = SAME_NONZERO_PATTERN;
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "MatDestroy_SuperLU"
273 PetscErrorCode MatDestroy_SuperLU(Mat A)
274 {
275   PetscErrorCode ierr;
276   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
277 
278   PetscFunctionBegin;
279   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
280     Destroy_SuperMatrix_Store(&lu->A);
281     Destroy_SuperMatrix_Store(&lu->B);
282     Destroy_SuperMatrix_Store(&lu->X);
283 
284     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
285     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
286     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
287     ierr = PetscFree(lu->R);CHKERRQ(ierr);
288     ierr = PetscFree(lu->C);CHKERRQ(ierr);
289     if ( lu->lwork >= 0 ) {
290       Destroy_SuperNode_Matrix(&lu->L);
291       Destroy_CompCol_Matrix(&lu->U);
292     }
293   }
294   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
295   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "MatView_SuperLU"
301 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
302 {
303   PetscErrorCode    ierr;
304   PetscTruth        iascii;
305   PetscViewerFormat format;
306   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);
307 
308   PetscFunctionBegin;
309   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
310 
311   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
312   if (iascii) {
313     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
314     if (format == PETSC_VIEWER_ASCII_INFO) {
315       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
316     }
317   }
318   PetscFunctionReturn(0);
319 }
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "MatAssemblyEnd_SuperLU"
323 PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
324   PetscErrorCode ierr;
325   Mat_SuperLU    *lu=(Mat_SuperLU*)(A->spptr);
326 
327   PetscFunctionBegin;
328   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
329   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
330   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
331   PetscFunctionReturn(0);
332 }
333 
334 
335 #undef __FUNCT__
336 #define __FUNCT__ "MatSolve_SuperLU_Private"
337 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
338 {
339   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
340   PetscScalar    *barray,*xarray;
341   PetscErrorCode ierr;
342   PetscInt       info,i;
343   SuperLUStat_t  stat;
344   PetscReal      ferr,berr;
345 
346   PetscFunctionBegin;
347   if ( lu->lwork == -1 ) {
348     PetscFunctionReturn(0);
349   }
350   lu->B.ncol = 1;   /* Set the number of right-hand side */
351   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
352   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
353 
354 #if defined(PETSC_USE_COMPLEX)
355   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
356   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
357 #else
358   ((DNformat*)lu->B.Store)->nzval = barray;
359   ((DNformat*)lu->X.Store)->nzval = xarray;
360 #endif
361 
362   /* Initialize the statistics variables. */
363   StatInit(&stat);
364 
365   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
366 #if defined(PETSC_USE_COMPLEX)
367   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
368            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
369            &lu->mem_usage, &stat, &info);
370 #else
371   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
372            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
373            &lu->mem_usage, &stat, &info);
374 #endif
375   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
376   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
377 
378   if ( !info || info == lu->A.ncol+1 ) {
379     if ( lu->options.IterRefine ) {
380       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
381       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
382       for (i = 0; i < 1; ++i)
383         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
384     }
385   } else if ( info > 0 ){
386     if ( lu->lwork == -1 ) {
387       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
388     } else {
389       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
390     }
391   } else if (info < 0){
392     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
393   }
394 
395   if ( lu->options.PrintStat ) {
396     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
397     StatPrint(&stat);
398   }
399   StatFree(&stat);
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "MatSolve_SuperLU"
405 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
406 {
407   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
408   PetscErrorCode ierr;
409 
410   PetscFunctionBegin;
411   lu->options.Trans = TRANS;
412   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "MatSolveTranspose_SuperLU"
418 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
419 {
420   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
421   PetscErrorCode ierr;
422 
423   PetscFunctionBegin;
424   lu->options.Trans = NOTRANS;
425   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 
430 /*
431    Note the r permutation is ignored
432 */
433 #undef __FUNCT__
434 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
435 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
436 {
437   Mat            B;
438   Mat_SuperLU    *lu;
439   PetscErrorCode ierr;
440   PetscInt       m=A->rmap.n,n=A->cmap.n,indx;
441   PetscTruth     flg;
442   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
443   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
444   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
445 
446   PetscFunctionBegin;
447   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
448   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
449   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
450   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
451 
452   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
453   B->ops->solve           = MatSolve_SuperLU;
454   B->ops->solvetranspose  = MatSolveTranspose_SuperLU;
455   B->factor               = FACTOR_LU;
456   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
457 
458   lu = (Mat_SuperLU*)(B->spptr);
459 
460   /* Set SuperLU options */
461     /* the default values for options argument:
462 	options.Fact = DOFACT;
463         options.Equil = YES;
464     	options.ColPerm = COLAMD;
465 	options.DiagPivotThresh = 1.0;
466     	options.Trans = NOTRANS;
467     	options.IterRefine = NOREFINE;
468     	options.SymmetricMode = NO;
469     	options.PivotGrowth = NO;
470     	options.ConditionNumber = NO;
471     	options.PrintStat = YES;
472     */
473   set_default_options(&lu->options);
474   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
475   lu->options.Equil = NO;
476   lu->options.PrintStat = NO;
477   lu->lwork = 0;   /* allocate space internally by system malloc */
478 
479   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
480   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
481   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
482   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
483   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
484   ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
485   if (flg) lu->options.SymmetricMode = YES;
486   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
487   ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
488   if (flg) lu->options.PivotGrowth = YES;
489   ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
490   if (flg) lu->options.ConditionNumber = YES;
491   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
492   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
493   ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
494   if (flg) lu->options.ReplaceTinyPivot = YES;
495   ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
496   if (flg) lu->options.PrintStat = YES;
497   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
498   if (lu->lwork > 0 ){
499     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
500   } else if (lu->lwork != 0 && lu->lwork != -1){
501     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
502     lu->lwork = 0;
503   }
504   PetscOptionsEnd();
505 
506 #ifdef SUPERLU2
507   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
508                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
509 #endif
510 
511   /* Allocate spaces (notice sizes are for the transpose) */
512   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
513   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
514   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
515   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr);
516   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr);
517 
518   /* create rhs and solution x without allocate space for .Store */
519 #if defined(PETSC_USE_COMPLEX)
520   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
521   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
522 #else
523   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
524   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
525 #endif
526 
527   lu->flg            = DIFFERENT_NONZERO_PATTERN;
528   lu->CleanUpSuperLU = PETSC_TRUE;
529 
530   *F = B;
531   ierr = PetscLogObjectMemory(B,(A->rmap.n+A->cmap.n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr);
532   PetscFunctionReturn(0);
533 }
534 
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "MatDuplicate_SuperLU"
538 PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
539   PetscErrorCode ierr;
540   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
541 
542   PetscFunctionBegin;
543   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
544   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
545   PetscFunctionReturn(0);
546 }
547 
548 
549 /*MC
550   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
551   via the external package SuperLU.
552 
553   If SuperLU is installed (see the manual for
554   instructions on how to declare the existence of external packages),
555   a matrix type can be constructed which invokes SuperLU solvers.
556   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
557 
558   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation() is
559   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
560   the MATSEQAIJ type without data copy.
561 
562   Options Database Keys:
563 + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
564 . -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
565                                     1: MMD applied to A'*A,
566                                     2: MMD applied to A'+A,
567                                     3: COLAMD, approximate minimum degree column ordering
568 . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
569                           choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
570 - -mat_superlu_printstat - print SuperLU statistics about the factorization
571 
572    Level: beginner
573 
574 .seealso: PCLU
575 M*/
576 
577 /*
578     Constructor for the new derived matrix class. It simply creates the base
579    matrix class and then adds the additional information/methods needed by SuperLU.
580 */
581 EXTERN_C_BEGIN
582 #undef __FUNCT__
583 #define __FUNCT__ "MatCreate_SuperLU"
584 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A)
585 {
586   PetscErrorCode ierr;
587 
588   PetscFunctionBegin;
589   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
590   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 EXTERN_C_END
594