xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1 
2 /*
3         Provides an interface to the SuperLU 3.0 sparse solver
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h"
7 
8 EXTERN_C_BEGIN
9 #if defined(PETSC_USE_COMPLEX)
10 #include "zsp_defs.h"
11 #else
12 #include "dsp_defs.h"
13 #endif
14 #include "util.h"
15 EXTERN_C_END
16 
17 typedef struct {
18   SuperMatrix       A,L,U,B,X;
19   superlu_options_t options;
20   int               *perm_c; /* column permutation vector */
21   int               *perm_r; /* row permutations from partial pivoting */
22   int               *etree;
23   double            *R, *C;
24   char              equed[1];
25   int               lwork;
26   void              *work;
27   double            rpg, rcond;
28   mem_usage_t       mem_usage;
29   MatStructure      flg;
30 
31   /* A few function pointers for inheritance */
32   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
33   PetscErrorCode (*MatView)(Mat,PetscViewer);
34   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
35   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
36   PetscErrorCode (*MatDestroy)(Mat);
37 
38   /* Flag to clean up (non-global) SuperLU objects during Destroy */
39   PetscTruth CleanUpSuperLU;
40 } Mat_SuperLU;
41 
42 
43 EXTERN PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
44 EXTERN PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
45 
46 EXTERN_C_BEGIN
47 EXTERN PetscErrorCode MatConvert_SuperLU_SeqAIJ(Mat,const MatType,Mat*);
48 EXTERN PetscErrorCode MatConvert_SeqAIJ_SuperLU(Mat,const MatType,Mat*);
49 EXTERN_C_END
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "MatDestroy_SuperLU"
53 PetscErrorCode MatDestroy_SuperLU(Mat A)
54 {
55   PetscErrorCode ierr;
56   Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
57 
58   PetscFunctionBegin;
59   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
60     Destroy_SuperMatrix_Store(&lu->A);
61     Destroy_SuperMatrix_Store(&lu->B);
62     Destroy_SuperMatrix_Store(&lu->X);
63 
64     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
65     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
66     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
67     ierr = PetscFree(lu->R);CHKERRQ(ierr);
68     ierr = PetscFree(lu->C);CHKERRQ(ierr);
69     if ( lu->lwork >= 0 ) {
70       Destroy_SuperNode_Matrix(&lu->L);
71       Destroy_CompCol_Matrix(&lu->U);
72     }
73   }
74   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
75   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "MatView_SuperLU"
81 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
82 {
83   PetscErrorCode ierr;
84   PetscTruth        iascii;
85   PetscViewerFormat format;
86   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);
87 
88   PetscFunctionBegin;
89   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
90 
91   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
92   if (iascii) {
93     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
94     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
95       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
96     }
97   }
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "MatAssemblyEnd_SuperLU"
103 PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
104   PetscErrorCode ierr;
105   Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr);
106 
107   PetscFunctionBegin;
108   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
109 
110   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
111   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
112   PetscFunctionReturn(0);
113 }
114 
115 /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */
116 #ifdef SuperLU2
117 #include "src/mat/impls/dense/seq/dense.h"
118 #undef __FUNCT__
119 #define __FUNCT__ "MatCreateNull_SuperLU"
120 PetscErrorCode MatCreateNull_SuperLU(Mat A,Mat *nullMat)
121 {
122   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
123   int           numRows = A->m,numCols = A->n;
124   SCformat      *Lstore;
125   int           numNullCols,size;
126   SuperLUStat_t stat;
127 #if defined(PETSC_USE_COMPLEX)
128   doublecomplex *nullVals,*workVals;
129 #else
130   PetscScalar   *nullVals,*workVals;
131 #endif
132   int           row,newRow,col,newCol,block,b;
133   PetscErrorCode 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) {
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 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
197 {
198   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
199   PetscScalar   *barray,*xarray;
200   PetscErrorCode ierr;
201   int           info,i;
202   SuperLUStat_t stat;
203   double        ferr,berr;
204 
205   PetscFunctionBegin;
206   if ( lu->lwork == -1 ) {
207     PetscFunctionReturn(0);
208   }
209   lu->B.ncol = 1;   /* Set the number of right-hand side */
210   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
211   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
212 
213 #if defined(PETSC_USE_COMPLEX)
214   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
215   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
216 #else
217   ((DNformat*)lu->B.Store)->nzval = barray;
218   ((DNformat*)lu->X.Store)->nzval = xarray;
219 #endif
220 
221   /* Initialize the statistics variables. */
222   StatInit(&stat);
223 
224   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
225   lu->options.Trans = TRANS;
226 #if defined(PETSC_USE_COMPLEX)
227   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
228            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
229            &lu->mem_usage, &stat, &info);
230 #else
231   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
232            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
233            &lu->mem_usage, &stat, &info);
234 #endif
235   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
236   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
237 
238   if ( !info || info == lu->A.ncol+1 ) {
239     if ( lu->options.IterRefine ) {
240       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
241       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
242       for (i = 0; i < 1; ++i)
243         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
244     }
245   } else if ( info > 0 ){
246     if ( lu->lwork == -1 ) {
247       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
248     } else {
249       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
250     }
251   } else if (info < 0){
252     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
253   }
254 
255   if ( lu->options.PrintStat ) {
256     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
257     StatPrint(&stat);
258   }
259   StatFree(&stat);
260   PetscFunctionReturn(0);
261 }
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
265 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
266 {
267   Mat_SeqAIJ    *aa = (Mat_SeqAIJ*)(A)->data;
268   Mat_SuperLU   *lu = (Mat_SuperLU*)(*F)->spptr;
269   PetscErrorCode ierr;
270   int           info;
271   SuperLUStat_t stat;
272   double        ferr, berr;
273   NCformat      *Ustore;
274   SCformat      *Lstore;
275 
276   PetscFunctionBegin;
277   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
278     lu->options.Fact = SamePattern;
279     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
280     Destroy_SuperMatrix_Store(&lu->A);
281     if ( lu->lwork >= 0 ) {
282       Destroy_SuperNode_Matrix(&lu->L);
283       Destroy_CompCol_Matrix(&lu->U);
284       lu->options.Fact = SamePattern;
285     }
286   }
287 
288   /* Create the SuperMatrix for lu->A=A^T:
289        Since SuperLU likes column-oriented matrices,we pass it the transpose,
290        and then solve A^T X = B in MatSolve(). */
291 #if defined(PETSC_USE_COMPLEX)
292   zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
293                            SLU_NC,SLU_Z,SLU_GE);
294 #else
295   dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i,
296                            SLU_NC,SLU_D,SLU_GE);
297 #endif
298 
299   /* Initialize the statistics variables. */
300   StatInit(&stat);
301 
302   /* Numerical factorization */
303   lu->B.ncol = 0;  /* Indicate not to solve the system */
304 #if defined(PETSC_USE_COMPLEX)
305    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
306            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
307            &lu->mem_usage, &stat, &info);
308 #else
309   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
310            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
311            &lu->mem_usage, &stat, &info);
312 #endif
313   if ( !info || info == lu->A.ncol+1 ) {
314     if ( lu->options.PivotGrowth )
315       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
316     if ( lu->options.ConditionNumber )
317       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
318   } else if ( info > 0 ){
319     if ( lu->lwork == -1 ) {
320       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
321     } else {
322       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
323     }
324   } else { /* info < 0 */
325     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
326   }
327 
328   if ( lu->options.PrintStat ) {
329     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
330     StatPrint(&stat);
331     Lstore = (SCformat *) lu->L.Store;
332     Ustore = (NCformat *) lu->U.Store;
333     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %d\n", Lstore->nnz);
334     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %d\n", Ustore->nnz);
335     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
336     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
337 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
338 	       lu->mem_usage.expansions);
339   }
340   StatFree(&stat);
341 
342   lu->flg = SAME_NONZERO_PATTERN;
343   PetscFunctionReturn(0);
344 }
345 
346 /*
347    Note the r permutation is ignored
348 */
349 #undef __FUNCT__
350 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
351 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
352 {
353   Mat          B;
354   Mat_SuperLU  *lu;
355   PetscErrorCode ierr;
356   int          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 PetscErrorCode 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 PetscErrorCode 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 = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
493   PetscFunctionReturn(0);
494 }
495 
496 EXTERN_C_BEGIN
497 #undef __FUNCT__
498 #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
499 PetscErrorCode MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
500   /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
501   /* to its base PETSc type, so we will ignore 'MatType type'. */
502   int                  ierr;
503   Mat                  B=*newmat;
504   Mat_SuperLU   *lu=(Mat_SuperLU *)A->spptr;
505 
506   PetscFunctionBegin;
507   if (B != A) {
508     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
509   }
510   /* Reset the original function pointers */
511   B->ops->duplicate        = lu->MatDuplicate;
512   B->ops->view             = lu->MatView;
513   B->ops->assemblyend      = lu->MatAssemblyEnd;
514   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
515   B->ops->destroy          = lu->MatDestroy;
516   /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
517   ierr = PetscFree(lu);CHKERRQ(ierr);
518   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
519   *newmat = B;
520   PetscFunctionReturn(0);
521 }
522 EXTERN_C_END
523 
524 EXTERN_C_BEGIN
525 #undef __FUNCT__
526 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
527 PetscErrorCode MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,Mat *newmat) {
528   /* This routine is only called to convert to MATSUPERLU */
529   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
530   int         ierr;
531   Mat         B=*newmat;
532   Mat_SuperLU *lu;
533 
534   PetscFunctionBegin;
535   if (B != A) {
536     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
537   }
538 
539   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
540   lu->MatDuplicate         = A->ops->duplicate;
541   lu->MatView              = A->ops->view;
542   lu->MatAssemblyEnd       = A->ops->assemblyend;
543   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
544   lu->MatDestroy           = A->ops->destroy;
545   lu->CleanUpSuperLU       = PETSC_FALSE;
546 
547   B->spptr                 = (void*)lu;
548   B->ops->duplicate        = MatDuplicate_SuperLU;
549   B->ops->view             = MatView_SuperLU;
550   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
551   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
552   B->ops->choleskyfactorsymbolic = 0;
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 PetscErrorCode MatCreate_SuperLU(Mat A)
596 {
597   PetscErrorCode ierr;
598 
599   PetscFunctionBegin;
600   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
601   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
602   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
603   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
604   PetscFunctionReturn(0);
605 }
606 EXTERN_C_END
607