xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 72ce74bd30d5d09d1371a3eb387704f4af395d1a)
1 
2 /*
3         Provides an interface to the UMFPACKv4.3 sparse solver
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h"
7 
8 EXTERN_C_BEGIN
9 #include "umfpack.h"
10 EXTERN_C_END
11 
12 typedef struct {
13   void         *Symbolic, *Numeric;
14   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
15   int          *Wi,*ai,*aj,*perm_c;
16   PetscScalar  *av;
17   MatStructure flg;
18   PetscTruth   PetscMatOdering;
19 
20   /* A few function pointers for inheritance */
21   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
22   int (*MatView)(Mat,PetscViewer);
23   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
24   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
25   int (*MatDestroy)(Mat);
26 
27   /* Flag to clean up UMFPACK objects during Destroy */
28   PetscTruth CleanUpUMFPACK;
29 } Mat_UMFPACK;
30 
31 EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
32 
33 EXTERN_C_BEGIN
34 #undef __FUNCT__
35 #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
36 PetscErrorCode MatConvert_UMFPACK_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
37   /* This routine is only called to convert an unfactored PETSc-UMFPACK matrix */
38   /* to its base PETSc type, so we will ignore 'MatType type'. */
39   PetscErrorCode ierr;
40   Mat         B=*newmat;
41   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
42 
43   PetscFunctionBegin;
44   if (B != A) {
45     /* This routine was inherited from SeqAIJ. */
46     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
47   }
48   /* Reset the original function pointers */
49   A->ops->duplicate        = lu->MatDuplicate;
50   A->ops->view             = lu->MatView;
51   A->ops->assemblyend      = lu->MatAssemblyEnd;
52   A->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
53   A->ops->destroy          = lu->MatDestroy;
54 
55   ierr = PetscFree(lu);CHKERRQ(ierr);
56   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
57   *newmat = B;
58   PetscFunctionReturn(0);
59 }
60 EXTERN_C_END
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "MatDestroy_UMFPACK"
64 PetscErrorCode MatDestroy_UMFPACK(Mat A)
65 {
66   PetscErrorCode ierr;
67   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
68 
69   PetscFunctionBegin;
70   if (lu->CleanUpUMFPACK) {
71     umfpack_di_free_symbolic(&lu->Symbolic);
72     umfpack_di_free_numeric(&lu->Numeric);
73     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
74     ierr = PetscFree(lu->W);CHKERRQ(ierr);
75     if (lu->PetscMatOdering) {
76       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
77     }
78   }
79   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
80   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "MatSolve_UMFPACK"
86 PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x) {
87   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
88   PetscScalar *av=lu->av,*ba,*xa;
89   PetscErrorCode ierr;
90   int          *ai=lu->ai,*aj=lu->aj,status;
91 
92   PetscFunctionBegin;
93   /* solve Ax = b by umfpack_di_wsolve */
94   /* ----------------------------------*/
95   ierr = VecGetArray(b,&ba);
96   ierr = VecGetArray(x,&xa);
97 
98   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
99   umfpack_di_report_info(lu->Control, lu->Info);
100   if (status < 0){
101     umfpack_di_report_status(lu->Control, status);
102     SETERRQ(1,"umfpack_di_wsolve failed");
103   }
104 
105   ierr = VecRestoreArray(b,&ba);
106   ierr = VecRestoreArray(x,&xa);
107 
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
113 PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,Mat *F) {
114   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
115   int         *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr;
116   PetscScalar *av=lu->av;
117 
118   PetscFunctionBegin;
119   /* numeric factorization of A' */
120   /* ----------------------------*/
121   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
122       umfpack_di_free_numeric(&lu->Numeric);
123   }
124   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
125   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
126   /* report numeric factorization of A' when Control[PRL] > 3 */
127   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control);
128 
129   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
130     /* allocate working space to be used by Solve */
131     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
132     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
133   }
134   lu->flg = SAME_NONZERO_PATTERN;
135   lu->CleanUpUMFPACK = PETSC_TRUE;
136   PetscFunctionReturn(0);
137 }
138 
139 /*
140    Note the r permutation is ignored
141 */
142 #undef __FUNCT__
143 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
144 PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
145   Mat         B;
146   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
147   Mat_UMFPACK *lu;
148   PetscErrorCode ierr;
149   int          m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ra,idx;
150   PetscScalar *av=mat->a;
151   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
152               *scale[]={"NONE","SUM","MAX"};
153   PetscTruth  flg;
154 
155   PetscFunctionBegin;
156   /* Create the factorization matrix F */
157   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
158   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
159   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
160 
161   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
162   B->ops->solve           = MatSolve_UMFPACK;
163   B->factor               = FACTOR_LU;
164   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
165 
166   lu = (Mat_UMFPACK*)(B->spptr);
167 
168   /* initializations */
169   /* ------------------------------------------------*/
170   /* get the default control parameters */
171   umfpack_di_defaults (lu->Control);
172   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
173   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
174 
175   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
176   /* Control parameters used by reporting routiones */
177   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
178 
179   /* Control parameters for symbolic factorization */
180   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
181   if (flg) {
182     switch (idx){
183     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
184     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
185     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
186     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
187     }
188   }
189   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);
190   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);
191   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);
192   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);
193   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);
194   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
195   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
196 
197   /* Control parameters used by numeric factorization */
198   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);
199   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);
200   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
201   if (flg) {
202     switch (idx){
203     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
204     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
205     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
206     }
207   }
208   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);
209   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);
210   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
211 
212   /* Control parameters used by solve */
213   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
214 
215   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
216   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
217   if (lu->PetscMatOdering) {
218     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
219     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
220     ierr = PetscMemcpy(lu->perm_c,ra,A->m*sizeof(int));CHKERRQ(ierr);
221     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
222   }
223   PetscOptionsEnd();
224 
225   /* print the control parameters */
226   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
227 
228   /* symbolic factorization of A' */
229   /* ---------------------------------------------------------------------- */
230   if (lu->PetscMatOdering) { /* use Petsc row ordering */
231     status = umfpack_di_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
232   } else { /* use Umfpack col ordering */
233     status = umfpack_di_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
234   }
235   if (status < 0){
236     umfpack_di_report_info(lu->Control, lu->Info);
237     umfpack_di_report_status(lu->Control, status);
238     SETERRQ(1,"umfpack_di_symbolic failed");
239   }
240   /* report sumbolic factorization of A' when Control[PRL] > 3 */
241   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control);
242 
243   lu->flg = DIFFERENT_NONZERO_PATTERN;
244   lu->ai  = ai;
245   lu->aj  = aj;
246   lu->av  = av;
247   lu->CleanUpUMFPACK = PETSC_TRUE;
248   *F = B;
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
254 PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) {
255   PetscErrorCode ierr;
256   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
257 
258   PetscFunctionBegin;
259   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
260   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
261   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
262   PetscFunctionReturn(0);
263 }
264 
265 /* used by -ksp_view */
266 #undef __FUNCT__
267 #define __FUNCT__ "MatFactorInfo_UMFPACK"
268 PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) {
269   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
270   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   /* check if matrix is UMFPACK type */
274   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
275 
276   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
277   /* Control parameters used by reporting routiones */
278   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
279 
280   /* Control parameters used by symbolic factorization */
281   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
282   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
283   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
284   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
285   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
286   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
287   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
288   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
289 
290   /* Control parameters used by numeric factorization */
291   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
292   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
293   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
294   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
295   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
296 
297   /* Control parameters used by solve */
298   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
299 
300   /* mat ordering */
301   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
302 
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "MatView_UMFPACK"
308 PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer) {
309   PetscErrorCode ierr;
310   PetscTruth        iascii;
311   PetscViewerFormat format;
312   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
313 
314   PetscFunctionBegin;
315   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
316 
317   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
318   if (iascii) {
319     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
320     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
321       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
322     }
323   }
324   PetscFunctionReturn(0);
325 }
326 
327 EXTERN_C_BEGIN
328 #undef __FUNCT__
329 #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
330 PetscErrorCode MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat) {
331   /* This routine is only called to convert to MATUMFPACK */
332   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
333   PetscErrorCode ierr;
334   Mat         B=*newmat;
335   Mat_UMFPACK *lu;
336 
337   PetscFunctionBegin;
338   if (B != A) {
339     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
340   }
341 
342   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
343   lu->MatDuplicate         = A->ops->duplicate;
344   lu->MatView              = A->ops->view;
345   lu->MatAssemblyEnd       = A->ops->assemblyend;
346   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
347   lu->MatDestroy           = A->ops->destroy;
348   lu->CleanUpUMFPACK       = PETSC_FALSE;
349 
350   B->spptr                 = (void*)lu;
351   B->ops->duplicate        = MatDuplicate_UMFPACK;
352   B->ops->view             = MatView_UMFPACK;
353   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
354   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
355   B->ops->destroy          = MatDestroy_UMFPACK;
356 
357   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
358                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
359   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
360                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
361   PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.");
362   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
363   *newmat = B;
364   PetscFunctionReturn(0);
365 }
366 EXTERN_C_END
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "MatDuplicate_UMFPACK"
370 PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) {
371   PetscErrorCode ierr;
372   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
373 
374   PetscFunctionBegin;
375   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
376   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 /*MC
381   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
382   via the external package UMFPACK.
383 
384   If UMFPACK is installed (see the manual for
385   instructions on how to declare the existence of external packages),
386   a matrix type can be constructed which invokes UMFPACK solvers.
387   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
388   This matrix type is only supported for double precision real.
389 
390   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
391   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
392   the MATSEQAIJ type without data copy.
393 
394   Consult UMFPACK documentation for more information about the Control parameters
395   which correspond to the options database keys below.
396 
397   Options Database Keys:
398 + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
399 . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
400 . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
401 . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
402 . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
403 . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
404 - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
405 
406    Level: beginner
407 
408 .seealso: PCLU
409 M*/
410 
411 EXTERN_C_BEGIN
412 #undef __FUNCT__
413 #define __FUNCT__ "MatCreate_UMFPACK"
414 PetscErrorCode MatCreate_UMFPACK(Mat A)
415 {
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */
420   ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr);
421   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
422   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr);
423   PetscFunctionReturn(0);
424 }
425 EXTERN_C_END
426