xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 56335db2edba962a6421bc90d84be8b461340ff3)
1 #define PETSCMAT_DLL
2 
3 /*
4    Provides an interface to the UMFPACKv5.1 sparse solver
5 
6    This interface uses the "UF_long version" of the UMFPACK API
7    (*_dl_* and *_zl_* routines, instead of *_di_* and *_zi_* routines)
8    so that UMFPACK can address more than 2Gb of memory on 64 bit
9    machines.
10 
11    If sizeof(UF_long) == 32 the interface only allocates the memory
12    necessary for UMFPACK's working arrays (W, Wi and possibly
13    perm_c). If sizeof(UF_long) == 64, in addition to allocating the
14    working arrays, the interface also has to re-allocate the matrix
15    index arrays (ai and aj, which must be stored as UF_long).
16 
17    The interface is implemented for both real and complex
18    arithmetic. Complex numbers are assumed to be in packed format,
19    which requires UMFPACK >= 4.4.
20 
21    We thank Christophe Geuzaine <geuzaine@acm.caltech.edu> for upgrading this interface to the UMFPACKv5.1
22 */
23 
24 #include "src/mat/impls/aij/seq/aij.h"
25 
26 EXTERN_C_BEGIN
27 #include "umfpack.h"
28 EXTERN_C_END
29 
30 typedef struct {
31   void         *Symbolic, *Numeric;
32   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
33   UF_long      *Wi,*ai,*aj,*perm_c;
34   PetscScalar  *av;
35   MatStructure flg;
36   PetscTruth   PetscMatOdering;
37 
38   /* A few function pointers for inheritance */
39   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
40   PetscErrorCode (*MatView)(Mat,PetscViewer);
41   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
42   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
43   PetscErrorCode (*MatDestroy)(Mat);
44 
45   /* Flag to clean up UMFPACK objects during Destroy */
46   PetscTruth CleanUpUMFPACK;
47 } Mat_UMFPACK;
48 
49 EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
50 
51 EXTERN_C_BEGIN
52 #undef __FUNCT__
53 #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
54 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
55 {
56   PetscErrorCode ierr;
57   Mat         B=*newmat;
58   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
59 
60   PetscFunctionBegin;
61   if (reuse == MAT_INITIAL_MATRIX) {
62     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
63   }
64   /* Reset the original function pointers */
65   B->ops->duplicate        = lu->MatDuplicate;
66   B->ops->view             = lu->MatView;
67   B->ops->assemblyend      = lu->MatAssemblyEnd;
68   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
69   B->ops->destroy          = lu->MatDestroy;
70 
71   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_umfpack_C","",PETSC_NULL);CHKERRQ(ierr);
72   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_umfpack_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
73 
74   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
75   *newmat = B;
76   PetscFunctionReturn(0);
77 }
78 EXTERN_C_END
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "MatDestroy_UMFPACK"
82 PetscErrorCode MatDestroy_UMFPACK(Mat A)
83 {
84   PetscErrorCode ierr;
85   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
86 
87   PetscFunctionBegin;
88   if (lu->CleanUpUMFPACK) {
89 #if defined(PETSC_USE_COMPLEX)
90     umfpack_zl_free_symbolic(&lu->Symbolic);
91     umfpack_zl_free_numeric(&lu->Numeric);
92 #else
93     umfpack_dl_free_symbolic(&lu->Symbolic);
94     umfpack_dl_free_numeric(&lu->Numeric);
95 #endif
96     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
97     ierr = PetscFree(lu->W);CHKERRQ(ierr);
98     if(sizeof(UF_long) != sizeof(int)){
99       ierr = PetscFree(lu->ai);CHKERRQ(ierr);
100       ierr = PetscFree(lu->aj);CHKERRQ(ierr);
101     }
102     if (lu->PetscMatOdering) {
103       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
104     }
105   }
106   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
107   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "MatSolve_UMFPACK"
113 PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
114 {
115   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
116   PetscScalar *av=lu->av,*ba,*xa;
117   PetscErrorCode ierr;
118   UF_long     *ai=lu->ai,*aj=lu->aj,status;
119 
120   PetscFunctionBegin;
121   /* solve Ax = b by umfpack_*_wsolve */
122   /* ----------------------------------*/
123 
124 #if defined(PETSC_USE_COMPLEX)
125   ierr = VecConjugate(b);
126 #endif
127 
128   ierr = VecGetArray(b,&ba);
129   ierr = VecGetArray(x,&xa);
130 
131 #if defined(PETSC_USE_COMPLEX)
132   status = umfpack_zl_wsolve(UMFPACK_At,ai,aj,(double*)av,NULL,(double*)xa,NULL,(double*)ba,NULL,
133 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
134   umfpack_zl_report_info(lu->Control, lu->Info);
135   if (status < 0){
136     umfpack_zl_report_status(lu->Control, status);
137     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_wsolve failed");
138   }
139 #else
140   status = umfpack_dl_wsolve(UMFPACK_At,ai,aj,av,xa,ba,
141 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
142   umfpack_dl_report_info(lu->Control, lu->Info);
143   if (status < 0){
144     umfpack_dl_report_status(lu->Control, status);
145     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_wsolve failed");
146   }
147 #endif
148 
149   ierr = VecRestoreArray(b,&ba);
150   ierr = VecRestoreArray(x,&xa);
151 
152 #if defined(PETSC_USE_COMPLEX)
153   ierr = VecConjugate(b);
154   ierr = VecConjugate(x);
155 #endif
156 
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
162 PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
163 {
164   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
165   PetscErrorCode ierr;
166   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status;
167   PetscScalar *av=lu->av;
168 
169   PetscFunctionBegin;
170   /* numeric factorization of A' */
171   /* ----------------------------*/
172 
173 #if defined(PETSC_USE_COMPLEX)
174   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
175     umfpack_zl_free_numeric(&lu->Numeric);
176   }
177   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
178   if (status < 0) {
179     umfpack_zl_report_status(lu->Control, status);
180     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
181   }
182   /* report numeric factorization of A' when Control[PRL] > 3 */
183   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
184 #else
185   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
186     umfpack_dl_free_numeric(&lu->Numeric);
187   }
188   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
189   if (status < 0) {
190     umfpack_zl_report_status(lu->Control, status);
191     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
192   }
193   /* report numeric factorization of A' when Control[PRL] > 3 */
194   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
195 #endif
196 
197   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
198     /* allocate working space to be used by Solve */
199     ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr);
200 #if defined(PETSC_USE_COMPLEX)
201     ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
202 #else
203     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
204 #endif
205   }
206 
207   lu->flg = SAME_NONZERO_PATTERN;
208   lu->CleanUpUMFPACK = PETSC_TRUE;
209   PetscFunctionReturn(0);
210 }
211 
212 /*
213    Note the r permutation is ignored
214 */
215 #undef __FUNCT__
216 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
217 PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
218 {
219   Mat         B;
220   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
221   Mat_UMFPACK *lu;
222   PetscErrorCode ierr;
223   int         i,m=A->rmap.n,n=A->cmap.n,*ra,idx;
224   UF_long     status;
225 
226   PetscScalar *av=mat->a;
227   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
228               *scale[]={"NONE","SUM","MAX"};
229   PetscTruth  flg;
230 
231   PetscFunctionBegin;
232   /* Create the factorization matrix F */
233   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
234   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
235   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
236   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
237 
238   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
239   B->ops->solve           = MatSolve_UMFPACK;
240   B->factor               = FACTOR_LU;
241   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
242 
243   lu = (Mat_UMFPACK*)(B->spptr);
244 
245   /* initializations */
246   /* ------------------------------------------------*/
247   /* get the default control parameters */
248 #if defined(PETSC_USE_COMPLEX)
249   umfpack_zl_defaults(lu->Control);
250 #else
251   umfpack_dl_defaults(lu->Control);
252 #endif
253   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
254   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
255 
256   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
257   /* Control parameters used by reporting routiones */
258   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
259 
260   /* Control parameters for symbolic factorization */
261   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
262   if (flg) {
263     switch (idx){
264     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
265     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
266     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
267     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
268     }
269   }
270   ierr = PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],PETSC_NULL);CHKERRQ(ierr);
271   ierr = PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],PETSC_NULL);CHKERRQ(ierr);
272   ierr = PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],PETSC_NULL);CHKERRQ(ierr);
273   ierr = PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],PETSC_NULL);CHKERRQ(ierr);
274   ierr = PetscOptionsReal("-mat_umfpack_2by2_tolerance","Control[UMFPACK_2BY2_TOLERANCE]","None",lu->Control[UMFPACK_2BY2_TOLERANCE],&lu->Control[UMFPACK_2BY2_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
275   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
276   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
277 
278   /* Control parameters used by numeric factorization */
279   ierr = PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
280   ierr = PetscOptionsReal("-mat_umfpack_sym_pivot_tolerance","Control[UMFPACK_SYM_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],&lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
281   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
282   if (flg) {
283     switch (idx){
284     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
285     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
286     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
287     }
288   }
289   ierr = PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr);
290   ierr = PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr);
291   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
292 
293   /* Control parameters used by solve */
294   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
295 
296   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
297   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
298   if (lu->PetscMatOdering) {
299     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
300     ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr);
301     /* we cannot simply memcpy on 64 bit archs */
302     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
303     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
304   }
305   PetscOptionsEnd();
306 
307   if(sizeof(UF_long) != sizeof(int)){
308     /* we cannot directly use mat->i and mat->j on 64 bit archs */
309     ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr);
310     ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr);
311     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
312     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
313   }
314   else{
315     lu->ai = (UF_long*)mat->i;
316     lu->aj = (UF_long*)mat->j;
317   }
318 
319   /* print the control parameters */
320 #if defined(PETSC_USE_COMPLEX)
321   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
322 #else
323   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
324 #endif
325 
326   /* symbolic factorization of A' */
327   /* ---------------------------------------------------------------------- */
328 #if defined(PETSC_USE_COMPLEX)
329   if (lu->PetscMatOdering) { /* use Petsc row ordering */
330     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
331 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
332   } else { /* use Umfpack col ordering */
333     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
334 				 &lu->Symbolic,lu->Control,lu->Info);
335   }
336   if (status < 0){
337     umfpack_zl_report_info(lu->Control, lu->Info);
338     umfpack_zl_report_status(lu->Control, status);
339     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
340   }
341   /* report sumbolic factorization of A' when Control[PRL] > 3 */
342   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
343 #else
344   if (lu->PetscMatOdering) { /* use Petsc row ordering */
345     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
346 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
347   } else { /* use Umfpack col ordering */
348     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
349 				 &lu->Symbolic,lu->Control,lu->Info);
350   }
351   if (status < 0){
352     umfpack_dl_report_info(lu->Control, lu->Info);
353     umfpack_dl_report_status(lu->Control, status);
354     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
355   }
356   /* report sumbolic factorization of A' when Control[PRL] > 3 */
357   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
358 #endif
359 
360   lu->flg = DIFFERENT_NONZERO_PATTERN;
361   lu->av  = av;
362   lu->CleanUpUMFPACK = PETSC_TRUE;
363   *F = B;
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
369 PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
370 {
371   PetscErrorCode ierr;
372   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
373 
374   PetscFunctionBegin;
375   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
376   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
377   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
378   PetscFunctionReturn(0);
379 }
380 
381 /* used by -ksp_view */
382 #undef __FUNCT__
383 #define __FUNCT__ "MatFactorInfo_UMFPACK"
384 PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
385 {
386   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
387   PetscErrorCode ierr;
388 
389   PetscFunctionBegin;
390   /* check if matrix is UMFPACK type */
391   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
392 
393   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
394   /* Control parameters used by reporting routiones */
395   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
396 
397   /* Control parameters used by symbolic factorization */
398   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
399   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
400   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
401   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
402   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
403   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
404   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
405   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
406 
407   /* Control parameters used by numeric factorization */
408   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
409   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
410   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
411   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
412   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
413 
414   /* Control parameters used by solve */
415   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
416 
417   /* mat ordering */
418   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
419 
420   PetscFunctionReturn(0);
421 }
422 
423 #undef __FUNCT__
424 #define __FUNCT__ "MatView_UMFPACK"
425 PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
426 {
427   PetscErrorCode ierr;
428   PetscTruth        iascii;
429   PetscViewerFormat format;
430   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
431 
432   PetscFunctionBegin;
433   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
434 
435   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
436   if (iascii) {
437     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
438     if (format == PETSC_VIEWER_ASCII_INFO) {
439       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
440     }
441   }
442   PetscFunctionReturn(0);
443 }
444 
445 EXTERN_C_BEGIN
446 #undef __FUNCT__
447 #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
448 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
449 {
450   PetscErrorCode ierr;
451   Mat            B=*newmat;
452   Mat_UMFPACK    *lu;
453 
454   PetscFunctionBegin;
455   if (reuse == MAT_INITIAL_MATRIX) {
456     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
457   }
458   B->ops->matsolve = 0;
459 
460   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
461   lu->MatDuplicate         = A->ops->duplicate;
462   lu->MatView              = A->ops->view;
463   lu->MatAssemblyEnd       = A->ops->assemblyend;
464   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
465   lu->MatDestroy           = A->ops->destroy;
466   lu->CleanUpUMFPACK       = PETSC_FALSE;
467 
468   B->spptr                 = (void*)lu;
469   B->ops->duplicate        = MatDuplicate_UMFPACK;
470   B->ops->view             = MatView_UMFPACK;
471   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
472   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
473   B->ops->destroy          = MatDestroy_UMFPACK;
474 
475   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
476                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
477   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
478                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
479   ierr = PetscInfo(A,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
480   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
481   *newmat = B;
482   PetscFunctionReturn(0);
483 }
484 EXTERN_C_END
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "MatDuplicate_UMFPACK"
488 PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
489 {
490   PetscErrorCode ierr;
491   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
492 
493   PetscFunctionBegin;
494   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
495   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
496   PetscFunctionReturn(0);
497 }
498 
499 /*MC
500   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
501   via the external package UMFPACK.
502 
503   If UMFPACK is installed (see the manual for
504   instructions on how to declare the existence of external packages),
505   a matrix type can be constructed which invokes UMFPACK solvers.
506   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
507 
508   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
509   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
510   the MATSEQAIJ type without data copy.
511 
512   Consult UMFPACK documentation for more information about the Control parameters
513   which correspond to the options database keys below.
514 
515   Options Database Keys:
516 + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
517 . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
518 . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
519 . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
520 . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
521 . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
522 - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
523 
524    Level: beginner
525 
526 .seealso: PCLU
527 M*/
528 
529 EXTERN_C_BEGIN
530 #undef __FUNCT__
531 #define __FUNCT__ "MatCreate_UMFPACK"
532 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A)
533 {
534   PetscErrorCode ierr;
535 
536   PetscFunctionBegin;
537   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
538   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
539   PetscFunctionReturn(0);
540 }
541 EXTERN_C_END
542 
543 
544 
545