xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision b2da817d11b8f98ec4a15fcfa65f2fa8cfedc826)
1 /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2 
3 /*
4         Provides an interface to the SuperLU 3.0 sparse solver
5 */
6 
7 #include "src/mat/impls/aij/seq/aij.h"
8 
9 EXTERN_C_BEGIN
10 #if defined(PETSC_USE_COMPLEX)
11 #include "zsp_defs.h"
12 #else
13 #include "dsp_defs.h"
14 #endif
15 #include "util.h"
16 EXTERN_C_END
17 
18 typedef struct {
19   SuperMatrix       A,L,U,B,X;
20   superlu_options_t options;
21   int               *perm_c; /* column permutation vector */
22   int               *perm_r; /* row permutations from partial pivoting */
23   int               *etree;
24   double            *R, *C;
25   char              equed[1];
26   int               lwork;
27   void              *work;
28   double            rpg, rcond;
29   mem_usage_t       mem_usage;
30   MatStructure      flg;
31 
32   /* A few function pointers for inheritance */
33   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
34   int (*MatView)(Mat,PetscViewer);
35   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
36   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
37   int (*MatDestroy)(Mat);
38 
39   /* Flag to clean up (non-global) SuperLU objects during Destroy */
40   PetscTruth CleanUpSuperLU;
41 } Mat_SuperLU;
42 
43 
44 EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer);
45 EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
46 
47 EXTERN_C_BEGIN
48 EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,const MatType,Mat*);
49 EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,const MatType,Mat*);
50 EXTERN_C_END
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "MatDestroy_SuperLU"
54 int MatDestroy_SuperLU(Mat A)
55 {
56   int         ierr;
57   Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
58 
59   PetscFunctionBegin;
60   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
61     Destroy_SuperMatrix_Store(&lu->A);
62     Destroy_SuperMatrix_Store(&lu->B);
63     Destroy_SuperMatrix_Store(&lu->X);
64 
65     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
66     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
67     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
68     ierr = PetscFree(lu->R);CHKERRQ(ierr);
69     ierr = PetscFree(lu->C);CHKERRQ(ierr);
70     if ( lu->lwork >= 0 ) {
71       Destroy_SuperNode_Matrix(&lu->L);
72       Destroy_CompCol_Matrix(&lu->U);
73     }
74   }
75   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
76   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "MatView_SuperLU"
82 int MatView_SuperLU(Mat A,PetscViewer viewer)
83 {
84   int               ierr;
85   PetscTruth        isascii;
86   PetscViewerFormat format;
87   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);
88 
89   PetscFunctionBegin;
90   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
91 
92   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
93   if (isascii) {
94     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
95     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
96       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
97     }
98   }
99   PetscFunctionReturn(0);
100 }
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "MatAssemblyEnd_SuperLU"
104 int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
105   int         ierr;
106   Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr);
107 
108   PetscFunctionBegin;
109   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
110 
111   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
112   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
113   PetscFunctionReturn(0);
114 }
115 
116 /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */
117 #ifdef SuperLU2
118 #include "src/mat/impls/dense/seq/dense.h"
119 #undef __FUNCT__
120 #define __FUNCT__ "MatCreateNull_SuperLU"
121 int MatCreateNull_SuperLU(Mat A,Mat *nullMat)
122 {
123   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
124   int           numRows = A->m,numCols = A->n;
125   SCformat      *Lstore;
126   int           numNullCols,size;
127   SuperLUStat_t stat;
128 #if defined(PETSC_USE_COMPLEX)
129   doublecomplex *nullVals,*workVals;
130 #else
131   PetscScalar   *nullVals,*workVals;
132 #endif
133   int           row,newRow,col,newCol,block,b,ierr;
134 
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(A,MAT_COOKIE);
137   PetscValidPointer(nullMat);
138   if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
139   numNullCols = numCols - numRows;
140   if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems");
141   /* Create the null matrix */
142   ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr);
143   if (numNullCols == 0) {
144     ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145     ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146     PetscFunctionReturn(0);
147   }
148 #if defined(PETSC_USE_COMPLEX)
149   nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
150 #else
151   nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
152 #endif
153   /* Copy in the columns */
154   Lstore = (SCformat*)lu->L.Store;
155   for(block = 0; block <= Lstore->nsuper; block++) {
156     newRow = Lstore->sup_to_col[block];
157     size   = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
158     for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
159       newCol = Lstore->rowind[col];
160       if (newCol >= numRows) {
161         for(b = 0; b < size; b++)
162 #if defined(PETSC_USE_COMPLEX)
163           nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
164 #else
165           nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
166 #endif
167       }
168     }
169   }
170   /* Permute rhs to form P^T_c B */
171   ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr);
172   for(b = 0; b < numNullCols; b++) {
173     for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
174     for(row = 0; row < numRows; row++) nullVals[b*numRows+row]   = workVals[row];
175   }
176   /* Backward solve the upper triangle A x = b */
177   for(b = 0; b < numNullCols; b++) {
178 #if defined(PETSC_USE_COMPLEX)
179     sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
180 #else
181     sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
182 #endif
183     if (ierr < 0)
184       SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr);
185   }
186   ierr = PetscFree(workVals);CHKERRQ(ierr);
187 
188   ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
189   ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190   PetscFunctionReturn(0);
191 }
192 #endif
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "MatSolve_SuperLU"
196 int MatSolve_SuperLU(Mat A,Vec b,Vec x)
197 {
198   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
199   PetscScalar   *barray,*xarray;
200   int           ierr,info,i;
201   SuperLUStat_t stat;
202   double        ferr,berr;
203 
204   PetscFunctionBegin;
205   if ( lu->lwork == -1 ) {
206     PetscFunctionReturn(0);
207   }
208   lu->B.ncol = 1;   /* Set the number of right-hand side */
209   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
210   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
211 
212 #if defined(PETSC_USE_COMPLEX)
213   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
214   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
215 #else
216   ((DNformat*)lu->B.Store)->nzval = barray;
217   ((DNformat*)lu->X.Store)->nzval = xarray;
218 #endif
219 
220   /* Initialize the statistics variables. */
221   StatInit(&stat);
222 
223   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
224   lu->options.Trans = TRANS;
225 #if defined(PETSC_USE_COMPLEX)
226   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
227            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
228            &lu->mem_usage, &stat, &info);
229 #else
230   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
231            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
232            &lu->mem_usage, &stat, &info);
233 #endif
234   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
235   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
236 
237   if ( info == 0 || info == lu->A.ncol+1 ) {
238     if ( lu->options.IterRefine ) {
239       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
240       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
241       for (i = 0; i < 1; ++i)
242         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
243     }
244   } else if ( info > 0 ){
245     if ( lu->lwork == -1 ) {
246       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
247     } else {
248       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
249     }
250   } else if (info < 0){
251     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
252   }
253 
254   if ( lu->options.PrintStat ) {
255     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
256     StatPrint(&stat);
257   }
258   StatFree(&stat);
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
264 int MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
265 {
266   Mat_SeqAIJ    *aa = (Mat_SeqAIJ*)(A)->data;
267   Mat_SuperLU   *lu = (Mat_SuperLU*)(*F)->spptr;
268   int           ierr,info;
269   PetscTruth    flag;
270   SuperLUStat_t stat;
271   double        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, &info);
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, &info);
311 #endif
312   if ( info == 0 || info == 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 ( info > 0 ){
318     if ( lu->lwork == -1 ) {
319       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
320     } else {
321       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
322     }
323   } else { /* info < 0 */
324     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
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 int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
351 {
352   Mat          B;
353   Mat_SuperLU  *lu;
354   int          ierr,m=A->m,n=A->n,indx;
355   PetscTruth   flg;
356   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
357   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
358   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
359 
360   PetscFunctionBegin;
361 
362   ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr);
363   ierr = MatSetType(B,MATSUPERLU);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(int),&lu->etree);CHKERRQ(ierr);
430   ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr);
431   ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
432   ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr);
433   ierr = PetscMalloc(m*sizeof(int),&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   PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU));
449   PetscFunctionReturn(0);
450 }
451 
452 /* used by -ksp_view */
453 #undef __FUNCT__
454 #define __FUNCT__ "MatFactorInfo_SuperLU"
455 int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
456 {
457   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
458   int               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 int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
485   int         ierr;
486   Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr;
487 
488   PetscFunctionBegin;
489   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
490   ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr);
491   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 EXTERN_C_BEGIN
496 #undef __FUNCT__
497 #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
498 int MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
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   int                  ierr;
502   Mat                  B=*newmat;
503   Mat_SuperLU   *lu=(Mat_SuperLU *)A->spptr;
504 
505   PetscFunctionBegin;
506   if (B != A) {
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   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
518   *newmat = B;
519   PetscFunctionReturn(0);
520 }
521 EXTERN_C_END
522 
523 EXTERN_C_BEGIN
524 #undef __FUNCT__
525 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
526 int MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,Mat *newmat) {
527   /* This routine is only called to convert to MATSUPERLU */
528   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
529   int         ierr;
530   Mat         B=*newmat;
531   Mat_SuperLU *lu;
532 
533   PetscFunctionBegin;
534   if (B != A) {
535     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
536   }
537 
538   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
539   lu->MatDuplicate         = A->ops->duplicate;
540   lu->MatView              = A->ops->view;
541   lu->MatAssemblyEnd       = A->ops->assemblyend;
542   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
543   lu->MatDestroy           = A->ops->destroy;
544   lu->CleanUpSuperLU       = PETSC_FALSE;
545 
546   B->spptr                 = (void*)lu;
547   B->ops->duplicate        = MatDuplicate_SuperLU;
548   B->ops->view             = MatView_SuperLU;
549   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
550   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
551   B->ops->destroy          = MatDestroy_SuperLU;
552 
553   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
554                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
555   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
556                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
557   PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.");
558   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
559   *newmat = B;
560   PetscFunctionReturn(0);
561 }
562 EXTERN_C_END
563 
564 /*MC
565   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
566   via the external package SuperLU.
567 
568   If SuperLU is installed (see the manual for
569   instructions on how to declare the existence of external packages),
570   a matrix type can be constructed which invokes SuperLU solvers.
571   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
572   This matrix type is only supported for double precision real.
573 
574   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
575   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
576   the MATSEQAIJ type without data copy.
577 
578   Options Database Keys:
579 + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
580 - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
581                                     1: MMD applied to A'*A,
582                                     2: MMD applied to A'+A,
583                                     3: COLAMD, approximate minimum degree column ordering
584 
585    Level: beginner
586 
587 .seealso: PCLU
588 M*/
589 
590 EXTERN_C_BEGIN
591 #undef __FUNCT__
592 #define __FUNCT__ "MatCreate_SuperLU"
593 int MatCreate_SuperLU(Mat A) {
594   int ierr;
595 
596   PetscFunctionBegin;
597   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
598   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
599   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
600   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
601   PetscFunctionReturn(0);
602 }
603 EXTERN_C_END
604