xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 0e5e90baed20fcf66590e8d3abe9cacb3e2155d0)
1 #define PETSCMAT_DLL
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   PetscInt          *perm_c; /* column permutation vector */
22   PetscInt          *perm_r; /* row permutations from partial pivoting */
23   PetscInt          *etree;
24   PetscReal         *R, *C;
25   char              equed[1];
26   PetscInt          lwork;
27   void              *work;
28   PetscReal         rpg, rcond;
29   mem_usage_t       mem_usage;
30   MatStructure      flg;
31 
32   /* A few function pointers for inheritance */
33   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
34   PetscErrorCode (*MatView)(Mat,PetscViewer);
35   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
36   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
37   PetscErrorCode (*MatDestroy)(Mat);
38 
39   /* Flag to clean up (non-global) SuperLU objects during Destroy */
40   PetscTruth CleanUpSuperLU;
41 } Mat_SuperLU;
42 
43 
44 EXTERN PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
45 EXTERN PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
46 
47 EXTERN_C_BEGIN
48 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat,MatType,MatReuse,Mat*);
49 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat,MatType,MatReuse,Mat*);
50 EXTERN_C_END
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "MatDestroy_SuperLU"
54 PetscErrorCode MatDestroy_SuperLU(Mat A)
55 {
56   PetscErrorCode 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,MAT_REUSE_MATRIX,&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 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
83 {
84   PetscErrorCode    ierr;
85   PetscTruth        iascii;
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,&iascii);CHKERRQ(ierr);
93   if (iascii) {
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 PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
105   PetscErrorCode ierr;
106   Mat_SuperLU    *lu=(Mat_SuperLU*)(A->spptr);
107 
108   PetscFunctionBegin;
109   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
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   PetscInt      numRows = A->m,numCols = A->n;
124   SCformat      *Lstore;
125   PetscInt      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   PetscInt           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,nullMat);CHKERRQ(ierr);
141   ierr = MatSetSizes(*nullMat,numRows,numNullCols,numRows,numNullCols);CHKERRQ(ierr);
142   ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr);
143   ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr);
144   if (!numNullCols) {
145     ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146     ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
147     PetscFunctionReturn(0);
148   }
149 #if defined(PETSC_USE_COMPLEX)
150   nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
151 #else
152   nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
153 #endif
154   /* Copy in the columns */
155   Lstore = (SCformat*)lu->L.Store;
156   for(block = 0; block <= Lstore->nsuper; block++) {
157     newRow = Lstore->sup_to_col[block];
158     size   = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
159     for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
160       newCol = Lstore->rowind[col];
161       if (newCol >= numRows) {
162         for(b = 0; b < size; b++)
163 #if defined(PETSC_USE_COMPLEX)
164           nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
165 #else
166           nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
167 #endif
168       }
169     }
170   }
171   /* Permute rhs to form P^T_c B */
172   ierr = PetscMalloc(numRows*sizeof(PetscReal),&workVals);CHKERRQ(ierr);
173   for(b = 0; b < numNullCols; b++) {
174     for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
175     for(row = 0; row < numRows; row++) nullVals[b*numRows+row]   = workVals[row];
176   }
177   /* Backward solve the upper triangle A x = b */
178   for(b = 0; b < numNullCols; b++) {
179 #if defined(PETSC_USE_COMPLEX)
180     sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
181 #else
182     sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
183 #endif
184     if (ierr < 0)
185       SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %D was invalid",-ierr);
186   }
187   ierr = PetscFree(workVals);CHKERRQ(ierr);
188 
189   ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190   ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191   PetscFunctionReturn(0);
192 }
193 #endif
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "MatSolve_SuperLU"
197 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
198 {
199   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
200   PetscScalar    *barray,*xarray;
201   PetscErrorCode ierr;
202   PetscInt       info,i;
203   SuperLUStat_t  stat;
204   PetscReal      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 || 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(PETSC_ERR_LIB, "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 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F)
267 {
268   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)(A)->data;
269   Mat_SuperLU    *lu = (Mat_SuperLU*)(*F)->spptr;
270   PetscErrorCode ierr;
271   PetscInt       sinfo;
272   SuperLUStat_t  stat;
273   PetscReal      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, &sinfo);
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, &sinfo);
313 #endif
314   if ( !sinfo || sinfo == 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 ( sinfo > 0 ){
320     if ( lu->lwork == -1 ) {
321       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
322     } else {
323       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
324     }
325   } else { /* sinfo < 0 */
326     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
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 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
353 {
354   Mat            B;
355   Mat_SuperLU    *lu;
356   PetscErrorCode ierr;
357   PetscInt       m=A->m,n=A->n,indx;
358   PetscTruth     flg;
359   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
360   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
361   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
362 
363   PetscFunctionBegin;
364   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
365   ierr = MatSetSizes(B,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
366   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
367   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
368 
369   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
370   B->ops->solve           = MatSolve_SuperLU;
371   B->factor               = FACTOR_LU;
372   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
373 
374   lu = (Mat_SuperLU*)(B->spptr);
375 
376   /* Set SuperLU options */
377     /* the default values for options argument:
378 	options.Fact = DOFACT;
379         options.Equil = YES;
380     	options.ColPerm = COLAMD;
381 	options.DiagPivotThresh = 1.0;
382     	options.Trans = NOTRANS;
383     	options.IterRefine = NOREFINE;
384     	options.SymmetricMode = NO;
385     	options.PivotGrowth = NO;
386     	options.ConditionNumber = NO;
387     	options.PrintStat = YES;
388     */
389   set_default_options(&lu->options);
390   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
391   lu->options.Equil = NO;
392   lu->options.PrintStat = NO;
393   lu->lwork = 0;   /* allocate space internally by system malloc */
394 
395   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
396   /*
397   ierr = PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
398   if (flg) lu->options.Equil = YES; -- not supported by the interface !!!
399   */
400   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
401   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
402   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
403   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
404   ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
405   if (flg) lu->options.SymmetricMode = YES;
406   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
407   ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
408   if (flg) lu->options.PivotGrowth = YES;
409   ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
410   if (flg) lu->options.ConditionNumber = YES;
411   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
412   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
413   ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
414   if (flg) lu->options.ReplaceTinyPivot = YES;
415   ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
416   if (flg) lu->options.PrintStat = YES;
417   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
418   if (lu->lwork > 0 ){
419     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
420   } else if (lu->lwork != 0 && lu->lwork != -1){
421     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
422     lu->lwork = 0;
423   }
424   PetscOptionsEnd();
425 
426 #ifdef SUPERLU2
427   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
428                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
429 #endif
430 
431   /* Allocate spaces (notice sizes are for the transpose) */
432   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
433   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
434   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
435   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr);
436   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr);
437 
438   /* create rhs and solution x without allocate space for .Store */
439 #if defined(PETSC_USE_COMPLEX)
440   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
441   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
442 #else
443   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
444   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
445 #endif
446 
447   lu->flg            = DIFFERENT_NONZERO_PATTERN;
448   lu->CleanUpSuperLU = PETSC_TRUE;
449 
450   *F = B;
451   ierr = PetscLogObjectMemory(B,(A->m+A->n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 /* used by -ksp_view */
456 #undef __FUNCT__
457 #define __FUNCT__ "MatFactorInfo_SuperLU"
458 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
459 {
460   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
461   PetscErrorCode    ierr;
462   superlu_options_t options;
463 
464   PetscFunctionBegin;
465   /* check if matrix is superlu_dist type */
466   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
467 
468   options = lu->options;
469   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
470   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
471   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
472   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
473   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
474   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
475   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
476   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
477   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
478   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
479   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
480   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
481 
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "MatDuplicate_SuperLU"
487 PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
488   PetscErrorCode ierr;
489   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
490 
491   PetscFunctionBegin;
492   ierr = (*lu->MatDuplicate)(A,op,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 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
501 {
502   /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
503   /* to its base PETSc type, so we will ignore 'MatType type'. */
504   PetscErrorCode ierr;
505   Mat            B=*newmat;
506   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
507 
508   PetscFunctionBegin;
509   if (reuse == MAT_INITIAL_MATRIX) {
510     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
511   }
512   /* Reset the original function pointers */
513   B->ops->duplicate        = lu->MatDuplicate;
514   B->ops->view             = lu->MatView;
515   B->ops->assemblyend      = lu->MatAssemblyEnd;
516   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
517   B->ops->destroy          = lu->MatDestroy;
518   /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
519   ierr = PetscFree(lu);CHKERRQ(ierr);
520 
521   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr);
522   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
523 
524   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
525   *newmat = B;
526   PetscFunctionReturn(0);
527 }
528 EXTERN_C_END
529 
530 EXTERN_C_BEGIN
531 #undef __FUNCT__
532 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
533 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat)
534 {
535   /* This routine is only called to convert to MATSUPERLU */
536   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
537   PetscErrorCode ierr;
538   Mat            B=*newmat;
539   Mat_SuperLU    *lu;
540 
541   PetscFunctionBegin;
542   if (reuse == MAT_INITIAL_MATRIX) {
543     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
544   }
545 
546   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
547   lu->MatDuplicate         = A->ops->duplicate;
548   lu->MatView              = A->ops->view;
549   lu->MatAssemblyEnd       = A->ops->assemblyend;
550   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
551   lu->MatDestroy           = A->ops->destroy;
552   lu->CleanUpSuperLU       = PETSC_FALSE;
553 
554   B->spptr                 = (void*)lu;
555   B->ops->duplicate        = MatDuplicate_SuperLU;
556   B->ops->view             = MatView_SuperLU;
557   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
558   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
559   B->ops->choleskyfactorsymbolic = 0;
560   B->ops->destroy          = MatDestroy_SuperLU;
561 
562   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
563                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
564   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
565                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
566   ierr = PetscLogInfo((0,"MatConvert_SeqAIJ_SuperLU:Using SuperLU for SeqAIJ LU factorization and solves.\n"));CHKERRQ(ierr);
567   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
568   *newmat = B;
569   PetscFunctionReturn(0);
570 }
571 EXTERN_C_END
572 
573 /*MC
574   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
575   via the external package SuperLU.
576 
577   If SuperLU is installed (see the manual for
578   instructions on how to declare the existence of external packages),
579   a matrix type can be constructed which invokes SuperLU solvers.
580   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
581 
582   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
583   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
584   the MATSEQAIJ type without data copy.
585 
586   Options Database Keys:
587 + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
588 - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
589                                     1: MMD applied to A'*A,
590                                     2: MMD applied to A'+A,
591                                     3: COLAMD, approximate minimum degree column ordering
592 
593    Level: beginner
594 
595 .seealso: PCLU
596 M*/
597 
598 EXTERN_C_BEGIN
599 #undef __FUNCT__
600 #define __FUNCT__ "MatCreate_SuperLU"
601 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A)
602 {
603   PetscErrorCode ierr;
604 
605   PetscFunctionBegin;
606   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
607   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
608   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
609   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
610   PetscFunctionReturn(0);
611 }
612 EXTERN_C_END
613