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