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