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