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