xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 7b6bb2c608b6fc6714ef38fda02c2dbb91c82665)
1 
2 /*
3    Provides an interface to the UMFPACKv5.1 sparse solver
4 
5    When build with PETSC_USE_64BIT_INDICES this will use UF_Long as the
6    integer type in UMFPACK, otherwise it will use int. This means
7    all integers in this file as simply declared as PetscInt. Also it means
8    that UMFPACK UL_Long version MUST be built with 64 bit integers when used.
9 
10 */
11 #include <../src/mat/impls/aij/seq/aij.h>
12 
13 #if defined(PETSC_USE_64BIT_INDICES)
14 #if defined(PETSC_USE_COMPLEX)
15 #define umfpack_UMF_free_symbolic   umfpack_zl_free_symbolic
16 #define umfpack_UMF_free_numeric    umfpack_zl_free_numeric
17 #define umfpack_UMF_wsolve          umfpack_zl_wsolve
18 #define umfpack_UMF_numeric         umfpack_zl_numeric
19 #define umfpack_UMF_report_numeric  umfpack_zl_report_numeric
20 #define umfpack_UMF_report_control  umfpack_zl_report_control
21 #define umfpack_UMF_report_status   umfpack_zl_report_status
22 #define umfpack_UMF_report_info     umfpack_zl_report_info
23 #define umfpack_UMF_report_symbolic umfpack_zl_report_symbolic
24 #define umfpack_UMF_qsymbolic       umfpack_zl_qsymbolic
25 #define umfpack_UMF_symbolic        umfpack_zl_symbolic
26 #define umfpack_UMF_defaults        umfpack_zl_defaults
27 
28 #else
29 #define umfpack_UMF_free_symbolic   umfpack_dl_free_symbolic
30 #define umfpack_UMF_free_numeric    umfpack_dl_free_numeric
31 #define umfpack_UMF_wsolve          umfpack_dl_wsolve
32 #define umfpack_UMF_numeric         umfpack_dl_numeric
33 #define umfpack_UMF_report_numeric  umfpack_dl_report_numeric
34 #define umfpack_UMF_report_control  umfpack_dl_report_control
35 #define umfpack_UMF_report_status   umfpack_dl_report_status
36 #define umfpack_UMF_report_info     umfpack_dl_report_info
37 #define umfpack_UMF_report_symbolic umfpack_dl_report_symbolic
38 #define umfpack_UMF_qsymbolic       umfpack_dl_qsymbolic
39 #define umfpack_UMF_symbolic        umfpack_dl_symbolic
40 #define umfpack_UMF_defaults        umfpack_dl_defaults
41 #endif
42 
43 #else
44 #if defined(PETSC_USE_COMPLEX)
45 #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
46 #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
47 #define umfpack_UMF_wsolve          umfpack_zi_wsolve
48 #define umfpack_UMF_numeric         umfpack_zi_numeric
49 #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
50 #define umfpack_UMF_report_control  umfpack_zi_report_control
51 #define umfpack_UMF_report_status   umfpack_zi_report_status
52 #define umfpack_UMF_report_info     umfpack_zi_report_info
53 #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
54 #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
55 #define umfpack_UMF_symbolic        umfpack_zi_symbolic
56 #define umfpack_UMF_defaults        umfpack_zi_defaults
57 
58 #else
59 #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
60 #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
61 #define umfpack_UMF_wsolve          umfpack_di_wsolve
62 #define umfpack_UMF_numeric         umfpack_di_numeric
63 #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
64 #define umfpack_UMF_report_control  umfpack_di_report_control
65 #define umfpack_UMF_report_status   umfpack_di_report_status
66 #define umfpack_UMF_report_info     umfpack_di_report_info
67 #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
68 #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
69 #define umfpack_UMF_symbolic        umfpack_di_symbolic
70 #define umfpack_UMF_defaults        umfpack_di_defaults
71 #endif
72 #endif
73 
74 
75 #define UF_long long long
76 #define UF_long_max LONG_LONG_MAX
77 #define UF_long_id "%lld"
78 
79 EXTERN_C_BEGIN
80 #include <umfpack.h>
81 EXTERN_C_END
82 
83 static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};
84 
85 typedef struct {
86   void         *Symbolic, *Numeric;
87   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
88   PetscInt      *Wi,*ai,*aj,*perm_c;
89   PetscScalar  *av;
90   MatStructure flg;
91   PetscBool    PetscMatOrdering;
92 
93   /* Flag to clean up UMFPACK objects during Destroy */
94   PetscBool  CleanUpUMFPACK;
95 } Mat_UMFPACK;
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MatDestroy_UMFPACK"
99 static PetscErrorCode MatDestroy_UMFPACK(Mat A)
100 {
101   PetscErrorCode ierr;
102   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
103 
104   PetscFunctionBegin;
105   if (lu && lu->CleanUpUMFPACK) {
106     umfpack_UMF_free_symbolic(&lu->Symbolic);
107     umfpack_UMF_free_numeric(&lu->Numeric);
108     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
109     ierr = PetscFree(lu->W);CHKERRQ(ierr);
110     if (lu->PetscMatOrdering) {
111       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
112     }
113   }
114   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
115   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "MatSolve_UMFPACK_Private"
121 static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
122 {
123   Mat_UMFPACK    *lu = (Mat_UMFPACK*)A->spptr;
124   PetscScalar    *av=lu->av,*ba,*xa;
125   PetscErrorCode ierr;
126   PetscInt       *ai=lu->ai,*aj=lu->aj,status;
127 
128   PetscFunctionBegin;
129   /* solve Ax = b by umfpack_*_wsolve */
130   /* ----------------------------------*/
131 
132   ierr = VecGetArray(b,&ba);
133   ierr = VecGetArray(x,&xa);
134 #if defined(PETSC_USE_COMPLEX)
135   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,
136                               lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
137 #else
138   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
139 #endif
140   umfpack_UMF_report_info(lu->Control, lu->Info);
141   if (status < 0){
142     umfpack_UMF_report_status(lu->Control, status);
143     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
144   }
145 
146   ierr = VecRestoreArray(b,&ba);CHKERRQ(ierr);
147   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatSolve_UMFPACK"
153 static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
154 {
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
159   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "MatSolveTranspose_UMFPACK"
165 static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
166 {
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
171   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
177 static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
178 {
179   Mat_UMFPACK *lu=(Mat_UMFPACK*)(F)->spptr;
180   PetscErrorCode ierr;
181   PetscInt     *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status;
182   PetscScalar *av=lu->av;
183 
184   PetscFunctionBegin;
185   /* numeric factorization of A' */
186   /* ----------------------------*/
187 
188   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
189     umfpack_UMF_free_numeric(&lu->Numeric);
190   }
191 #if defined(PETSC_USE_COMPLEX)
192   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
193 #else
194   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
195 #endif
196   if (status < 0) {
197     umfpack_UMF_report_status(lu->Control, status);
198     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
199   }
200   /* report numeric factorization of A' when Control[PRL] > 3 */
201   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
202 
203   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
204     /* allocate working space to be used by Solve */
205     ierr = PetscMalloc(m * sizeof(PetscInt), &lu->Wi);CHKERRQ(ierr);
206     ierr = PetscMalloc(5*m * sizeof(PetscScalar), &lu->W);CHKERRQ(ierr);
207   }
208 
209   lu->flg = SAME_NONZERO_PATTERN;
210   lu->CleanUpUMFPACK = PETSC_TRUE;
211   F->ops->solve          = MatSolve_UMFPACK;
212   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
213   PetscFunctionReturn(0);
214 }
215 
216 /*
217    Note the r permutation is ignored
218 */
219 #undef __FUNCT__
220 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
221 static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
222 {
223   Mat_SeqAIJ     *mat=(Mat_SeqAIJ*)A->data;
224   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
225   PetscErrorCode ierr;
226   PetscInt       i,m=A->rmap->n,n=A->cmap->n;
227   const PetscInt *ra;
228   PetscInt        status;
229   PetscScalar    *av=mat->a;
230 
231   PetscFunctionBegin;
232   if (lu->PetscMatOrdering) {
233     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
234     ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
235     /* we cannot simply memcpy on 64 bit archs */
236     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
237     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
238   }
239 
240   lu->ai = mat->i;
241   lu->aj = mat->j;
242 
243   /* print the control parameters */
244   if(lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
245 
246   /* symbolic factorization of A' */
247   /* ---------------------------------------------------------------------- */
248   if (lu->PetscMatOrdering) { /* use Petsc row ordering */
249 #if !defined(PETSC_USE_COMPLEX)
250     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
251 #else
252     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,NULL,NULL,
253                                    lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
254 #endif
255   } else { /* use Umfpack col ordering */
256 #if !defined(PETSC_USE_COMPLEX)
257     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,av,&lu->Symbolic,lu->Control,lu->Info);
258 #else
259     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
260 #endif
261   }
262   if (status < 0){
263     umfpack_UMF_report_info(lu->Control, lu->Info);
264     umfpack_UMF_report_status(lu->Control, status);
265     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
266   }
267   /* report sumbolic factorization of A' when Control[PRL] > 3 */
268   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
269 
270   lu->flg = DIFFERENT_NONZERO_PATTERN;
271   lu->av  = av;
272   lu->CleanUpUMFPACK = PETSC_TRUE;
273   (F)->ops->lufactornumeric  = MatLUFactorNumeric_UMFPACK;
274   PetscFunctionReturn(0);
275 }
276 
277 #undef __FUNCT__
278 #define __FUNCT__ "MatFactorInfo_UMFPACK"
279 static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
280 {
281   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
282   PetscErrorCode ierr;
283 
284   PetscFunctionBegin;
285   /* check if matrix is UMFPACK type */
286   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
287 
288   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
289   /* Control parameters used by reporting routiones */
290   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
291 
292   /* Control parameters used by symbolic factorization */
293   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
294   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
295   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
296   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
297   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
298   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
299   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
300 
301   /* Control parameters used by numeric factorization */
302   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
303   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
304   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
305   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
306   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
307 
308   /* Control parameters used by solve */
309   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
310 
311   /* mat ordering */
312   if (!lu->PetscMatOrdering) {
313     ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr);
314   }
315 
316   PetscFunctionReturn(0);
317 }
318 
319 #undef __FUNCT__
320 #define __FUNCT__ "MatView_UMFPACK"
321 static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
322 {
323   PetscErrorCode    ierr;
324   PetscBool         iascii;
325   PetscViewerFormat format;
326 
327   PetscFunctionBegin;
328   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
329 
330   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
331   if (iascii) {
332     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
333     if (format == PETSC_VIEWER_ASCII_INFO) {
334       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
335     }
336   }
337   PetscFunctionReturn(0);
338 }
339 
340 EXTERN_C_BEGIN
341 #undef __FUNCT__
342 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
343 PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
344 {
345   PetscFunctionBegin;
346   *type = MATSOLVERUMFPACK;
347   PetscFunctionReturn(0);
348 }
349 EXTERN_C_END
350 
351 
352 /*MC
353   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
354   via the external package UMFPACK.
355 
356   ./configure --download-umfpack to install PETSc to use UMFPACK
357 
358   Consult UMFPACK documentation for more information about the Control parameters
359   which correspond to the options database keys below.
360 
361   Options Database Keys:
362 + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
363 . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
364 . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
365 . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW]
366 . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE]
367 . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
368 . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE]
369 . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ]
370 . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE]
371 . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
372 .  -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE]
373 .  -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX
374 . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
375 .  -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL]
376 - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
377 
378    Level: beginner
379 
380 .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
381 M*/
382 EXTERN_C_BEGIN
383 #undef __FUNCT__
384 #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
385 PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
386 {
387   Mat            B;
388   Mat_UMFPACK    *lu;
389   PetscErrorCode ierr;
390   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
391 
392   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"},
393                  *scale[]={"NONE","SUM","MAX"};
394   PetscBool      flg;
395 
396   PetscFunctionBegin;
397   /* Create the factorization matrix F */
398   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
399   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
400   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
401   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
402   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
403   B->spptr                 = lu;
404   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
405   B->ops->destroy          = MatDestroy_UMFPACK;
406   B->ops->view             = MatView_UMFPACK;
407   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_umfpack",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
408   B->factortype            = MAT_FACTOR_LU;
409   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
410   B->preallocated          = PETSC_TRUE;
411 
412   /* initializations */
413   /* ------------------------------------------------*/
414   /* get the default control parameters */
415   umfpack_UMF_defaults(lu->Control);
416   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
417   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
418 
419   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
420   /* Control parameters used by reporting routiones */
421   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
422 
423   /* Control parameters for symbolic factorization */
424   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr);
425   if (flg) {
426     switch (idx){
427     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
428     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
429     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
430     }
431   }
432   ierr = PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,sizeof UmfpackOrderingTypes/sizeof UmfpackOrderingTypes[0],UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg);CHKERRQ(ierr);
433   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
434   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);
435   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);
436   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);
437   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);
438   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
439   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
440 
441   /* Control parameters used by numeric factorization */
442   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);
443   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);
444   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
445   if (flg) {
446     switch (idx){
447     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
448     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
449     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
450     }
451   }
452   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);
453   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);
454   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
455 
456   /* Control parameters used by solve */
457   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
458 
459   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
460   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOrdering);CHKERRQ(ierr);
461   PetscOptionsEnd();
462   *F = B;
463   PetscFunctionReturn(0);
464 }
465 EXTERN_C_END
466 
467