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