xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 01bab6536bbd56ce9fcee44f433cebab9e3b5e5f)
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 #endif
155 
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
161 PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
162 {
163   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
164   PetscErrorCode ierr;
165   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status;
166   PetscScalar *av=lu->av;
167 
168   PetscFunctionBegin;
169   /* numeric factorization of A' */
170   /* ----------------------------*/
171 
172 #if defined(PETSC_USE_COMPLEX)
173   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
174     umfpack_zl_free_numeric(&lu->Numeric);
175   }
176   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
177   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
178   /* report numeric factorization of A' when Control[PRL] > 3 */
179   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
180 #else
181   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
182     umfpack_dl_free_numeric(&lu->Numeric);
183   }
184   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
185   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
186   /* report numeric factorization of A' when Control[PRL] > 3 */
187   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
188 #endif
189 
190   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
191     /* allocate working space to be used by Solve */
192     ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr);
193 #if defined(PETSC_USE_COMPLEX)
194     ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
195 #else
196     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
197 #endif
198   }
199 
200   lu->flg = SAME_NONZERO_PATTERN;
201   lu->CleanUpUMFPACK = PETSC_TRUE;
202   PetscFunctionReturn(0);
203 }
204 
205 /*
206    Note the r permutation is ignored
207 */
208 #undef __FUNCT__
209 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
210 PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
211 {
212   Mat         B;
213   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
214   Mat_UMFPACK *lu;
215   PetscErrorCode ierr;
216   int         i,m=A->rmap.n,n=A->cmap.n,*ra,idx;
217   UF_long     status;
218 
219   PetscScalar *av=mat->a;
220   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
221               *scale[]={"NONE","SUM","MAX"};
222   PetscTruth  flg;
223 
224   PetscFunctionBegin;
225   /* Create the factorization matrix F */
226   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
227   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
228   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
229   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
230 
231   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
232   B->ops->solve           = MatSolve_UMFPACK;
233   B->factor               = FACTOR_LU;
234   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
235 
236   lu = (Mat_UMFPACK*)(B->spptr);
237 
238   /* initializations */
239   /* ------------------------------------------------*/
240   /* get the default control parameters */
241 #if defined(PETSC_USE_COMPLEX)
242   umfpack_zl_defaults(lu->Control);
243 #else
244   umfpack_dl_defaults(lu->Control);
245 #endif
246   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
247   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
248 
249   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
250   /* Control parameters used by reporting routiones */
251   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
252 
253   /* Control parameters for symbolic factorization */
254   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
255   if (flg) {
256     switch (idx){
257     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
258     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
259     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
260     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
261     }
262   }
263   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);
264   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);
265   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);
266   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);
267   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);
268   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
269   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
270 
271   /* Control parameters used by numeric factorization */
272   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);
273   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);
274   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
275   if (flg) {
276     switch (idx){
277     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
278     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
279     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
280     }
281   }
282   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);
283   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);
284   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
285 
286   /* Control parameters used by solve */
287   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
288 
289   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
290   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
291   if (lu->PetscMatOdering) {
292     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
293     ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr);
294     /* we cannot simply memcpy on 64 bit archs */
295     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
296     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
297   }
298   PetscOptionsEnd();
299 
300   if(sizeof(UF_long) != sizeof(int)){
301     /* we cannot directly use mat->i and mat->j on 64 bit archs */
302     ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr);
303     ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr);
304     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
305     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
306   }
307   else{
308     lu->ai = (UF_long*)mat->i;
309     lu->aj = (UF_long*)mat->j;
310   }
311 
312   /* print the control parameters */
313 #if defined(PETSC_USE_COMPLEX)
314   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
315 #else
316   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
317 #endif
318 
319   /* symbolic factorization of A' */
320   /* ---------------------------------------------------------------------- */
321 #if defined(PETSC_USE_COMPLEX)
322   if (lu->PetscMatOdering) { /* use Petsc row ordering */
323     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
324 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
325   } else { /* use Umfpack col ordering */
326     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
327 				 &lu->Symbolic,lu->Control,lu->Info);
328   }
329   if (status < 0){
330     umfpack_zl_report_info(lu->Control, lu->Info);
331     umfpack_zl_report_status(lu->Control, status);
332     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
333   }
334   /* report sumbolic factorization of A' when Control[PRL] > 3 */
335   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
336 #else
337   if (lu->PetscMatOdering) { /* use Petsc row ordering */
338     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
339 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
340   } else { /* use Umfpack col ordering */
341     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
342 				 &lu->Symbolic,lu->Control,lu->Info);
343   }
344   if (status < 0){
345     umfpack_dl_report_info(lu->Control, lu->Info);
346     umfpack_dl_report_status(lu->Control, status);
347     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
348   }
349   /* report sumbolic factorization of A' when Control[PRL] > 3 */
350   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
351 #endif
352 
353   lu->flg = DIFFERENT_NONZERO_PATTERN;
354   lu->av  = av;
355   lu->CleanUpUMFPACK = PETSC_TRUE;
356   *F = B;
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
362 PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
363 {
364   PetscErrorCode ierr;
365   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
366 
367   PetscFunctionBegin;
368   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
369   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
370   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
371   PetscFunctionReturn(0);
372 }
373 
374 /* used by -ksp_view */
375 #undef __FUNCT__
376 #define __FUNCT__ "MatFactorInfo_UMFPACK"
377 PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
378 {
379   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   /* check if matrix is UMFPACK type */
384   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
385 
386   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
387   /* Control parameters used by reporting routiones */
388   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
389 
390   /* Control parameters used by symbolic factorization */
391   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
392   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
393   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
394   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
395   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
396   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
397   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
398   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
399 
400   /* Control parameters used by numeric factorization */
401   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
402   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
403   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
404   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
405   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
406 
407   /* Control parameters used by solve */
408   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
409 
410   /* mat ordering */
411   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
412 
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "MatView_UMFPACK"
418 PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
419 {
420   PetscErrorCode ierr;
421   PetscTruth        iascii;
422   PetscViewerFormat format;
423   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
424 
425   PetscFunctionBegin;
426   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
427 
428   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
429   if (iascii) {
430     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
431     if (format == PETSC_VIEWER_ASCII_INFO) {
432       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
433     }
434   }
435   PetscFunctionReturn(0);
436 }
437 
438 EXTERN_C_BEGIN
439 #undef __FUNCT__
440 #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
441 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
442 {
443   PetscErrorCode ierr;
444   Mat            B=*newmat;
445   Mat_UMFPACK    *lu;
446 
447   PetscFunctionBegin;
448   if (reuse == MAT_INITIAL_MATRIX) {
449     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
450   }
451 
452   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
453   lu->MatDuplicate         = A->ops->duplicate;
454   lu->MatView              = A->ops->view;
455   lu->MatAssemblyEnd       = A->ops->assemblyend;
456   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
457   lu->MatDestroy           = A->ops->destroy;
458   lu->CleanUpUMFPACK       = PETSC_FALSE;
459 
460   B->spptr                 = (void*)lu;
461   B->ops->duplicate        = MatDuplicate_UMFPACK;
462   B->ops->view             = MatView_UMFPACK;
463   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
464   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
465   B->ops->destroy          = MatDestroy_UMFPACK;
466 
467   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
468                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
469   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
470                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
471   ierr = PetscInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
472   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
473   *newmat = B;
474   PetscFunctionReturn(0);
475 }
476 EXTERN_C_END
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "MatDuplicate_UMFPACK"
480 PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
481 {
482   PetscErrorCode ierr;
483   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
484 
485   PetscFunctionBegin;
486   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
487   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
488   PetscFunctionReturn(0);
489 }
490 
491 /*MC
492   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
493   via the external package UMFPACK.
494 
495   If UMFPACK is installed (see the manual for
496   instructions on how to declare the existence of external packages),
497   a matrix type can be constructed which invokes UMFPACK solvers.
498   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
499 
500   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
501   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
502   the MATSEQAIJ type without data copy.
503 
504   Consult UMFPACK documentation for more information about the Control parameters
505   which correspond to the options database keys below.
506 
507   Options Database Keys:
508 + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
509 . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
510 . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
511 . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
512 . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
513 . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
514 - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
515 
516    Level: beginner
517 
518 .seealso: PCLU
519 M*/
520 
521 EXTERN_C_BEGIN
522 #undef __FUNCT__
523 #define __FUNCT__ "MatCreate_UMFPACK"
524 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A)
525 {
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
530   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
531   PetscFunctionReturn(0);
532 }
533 EXTERN_C_END
534 
535 
536 
537