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