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