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