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