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