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