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