xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 1677d0b88caff045fb399525ae6e8dad78c032c0)
1*1677d0b8SKris Buschelman /*$Id: umfpack.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2*1677d0b8SKris Buschelman 
3*1677d0b8SKris Buschelman /*
4*1677d0b8SKris Buschelman         Provides an interface to the UMFPACK sparse solver
5*1677d0b8SKris Buschelman */
6*1677d0b8SKris Buschelman 
7*1677d0b8SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
8*1677d0b8SKris Buschelman 
9*1677d0b8SKris Buschelman EXTERN_C_BEGIN
10*1677d0b8SKris Buschelman #include "umfpack.h"
11*1677d0b8SKris Buschelman EXTERN_C_END
12*1677d0b8SKris Buschelman 
13*1677d0b8SKris Buschelman typedef struct {
14*1677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
15*1677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
16*1677d0b8SKris Buschelman   int          *Wi,*ai,*aj,*perm_c;
17*1677d0b8SKris Buschelman   PetscScalar  *av;
18*1677d0b8SKris Buschelman   MatStructure flg;
19*1677d0b8SKris Buschelman   PetscTruth   PetscMatOdering;
20*1677d0b8SKris Buschelman 
21*1677d0b8SKris Buschelman   /* A few function pointers for inheritance */
22*1677d0b8SKris Buschelman   int (*MatView)(Mat,PetscViewer);
23*1677d0b8SKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
24*1677d0b8SKris Buschelman   int (*MatDestroy)(Mat);
25*1677d0b8SKris Buschelman 
26*1677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
27*1677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
28*1677d0b8SKris Buschelman } Mat_SeqAIJ_UMFPACK;
29*1677d0b8SKris Buschelman 
30*1677d0b8SKris Buschelman #undef __FUNCT__
31*1677d0b8SKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_UMFPACK"
32*1677d0b8SKris Buschelman int MatDestroy_SeqAIJ_UMFPACK(Mat A)
33*1677d0b8SKris Buschelman {
34*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)A->spptr;
35*1677d0b8SKris Buschelman   int                ierr,(*destroy)(Mat);
36*1677d0b8SKris Buschelman 
37*1677d0b8SKris Buschelman   PetscFunctionBegin;
38*1677d0b8SKris Buschelman   if (CleanUpUMFPACK) {
39*1677d0b8SKris Buschelman     umfpack_di_free_symbolic(&lu->Symbolic) ;
40*1677d0b8SKris Buschelman     umfpack_di_free_numeric(&lu->Numeric) ;
41*1677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
42*1677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
43*1677d0b8SKris Buschelman 
44*1677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
45*1677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
46*1677d0b8SKris Buschelman     }
47*1677d0b8SKris Buschelman   }
48*1677d0b8SKris Buschelman 
49*1677d0b8SKris Buschelman   destroy = lu->MatDestroy;
50*1677d0b8SKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
51*1677d0b8SKris Buschelman   ierr = (*destroy)(A);CHKERRQ(ierr);
52*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
53*1677d0b8SKris Buschelman }
54*1677d0b8SKris Buschelman 
55*1677d0b8SKris Buschelman #undef __FUNCT__
56*1677d0b8SKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_UMFPACK"
57*1677d0b8SKris Buschelman int MatView_SeqAIJ_UMFPACK(Mat A,PetscViewer viewer)
58*1677d0b8SKris Buschelman {
59*1677d0b8SKris Buschelman   int                 ierr;
60*1677d0b8SKris Buschelman   PetscTruth          isascii;
61*1677d0b8SKris Buschelman   PetscViewerFormat   format;
62*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK  *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr);
63*1677d0b8SKris Buschelman 
64*1677d0b8SKris Buschelman   PetscFunctionBegin;
65*1677d0b8SKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
66*1677d0b8SKris Buschelman 
67*1677d0b8SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
68*1677d0b8SKris Buschelman   if (isascii) {
69*1677d0b8SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
70*1677d0b8SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
71*1677d0b8SKris Buschelman       ierr = MatSeqAIJFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
72*1677d0b8SKris Buschelman     }
73*1677d0b8SKris Buschelman   }
74*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
75*1677d0b8SKris Buschelman }
76*1677d0b8SKris Buschelman 
77*1677d0b8SKris Buschelman #undef __FUNCT__
78*1677d0b8SKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_UMFPACK"
79*1677d0b8SKris Buschelman int MatAssemblyEnd_SeqAIJ_UMFPACK(Mat A,MatAssemblyType mode) {
80*1677d0b8SKris Buschelman   int                ierr;
81*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr);
82*1677d0b8SKris Buschelman 
83*1677d0b8SKris Buschelman   PetscFunctionBegin;
84*1677d0b8SKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
85*1677d0b8SKris Buschelman   ierr = MatUseUMFPACK_SeqAIJ(A);CHKERRQ(ierr);
86*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
87*1677d0b8SKris Buschelman }
88*1677d0b8SKris Buschelman 
89*1677d0b8SKris Buschelman #undef __FUNCT__
90*1677d0b8SKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_UMFPACK"
91*1677d0b8SKris Buschelman int MatSolve_SeqAIJ_UMFPACK(Mat A,Vec b,Vec x)
92*1677d0b8SKris Buschelman {
93*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)A->spptr;
94*1677d0b8SKris Buschelman   PetscScalar        *av=lu->av,*ba,*xa;
95*1677d0b8SKris Buschelman   int                ierr,*ai=lu->ai,*aj=lu->aj,status;
96*1677d0b8SKris Buschelman 
97*1677d0b8SKris Buschelman   PetscFunctionBegin;
98*1677d0b8SKris Buschelman   /* solve Ax = b by umfpack_di_wsolve */
99*1677d0b8SKris Buschelman   /* ----------------------------------*/
100*1677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
101*1677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
102*1677d0b8SKris Buschelman 
103*1677d0b8SKris Buschelman   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
104*1677d0b8SKris Buschelman   umfpack_di_report_info(lu->Control, lu->Info);
105*1677d0b8SKris Buschelman   if (status < 0){
106*1677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
107*1677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_wsolve failed") ;
108*1677d0b8SKris Buschelman   }
109*1677d0b8SKris Buschelman 
110*1677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
111*1677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
112*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
113*1677d0b8SKris Buschelman }
114*1677d0b8SKris Buschelman 
115*1677d0b8SKris Buschelman #undef __FUNCT__
116*1677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_UMFPACK"
117*1677d0b8SKris Buschelman int MatLUFactorNumeric_SeqAIJ_UMFPACK(Mat A,Mat *F)
118*1677d0b8SKris Buschelman {
119*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)(*F)->spptr;
120*1677d0b8SKris Buschelman   int                *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr;
121*1677d0b8SKris Buschelman   PetscScalar        *av=lu->av;
122*1677d0b8SKris Buschelman 
123*1677d0b8SKris Buschelman   PetscFunctionBegin;
124*1677d0b8SKris Buschelman   /* numeric factorization of A' */
125*1677d0b8SKris Buschelman   /* ----------------------------*/
126*1677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ;
127*1677d0b8SKris Buschelman   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
128*1677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
129*1677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ;
130*1677d0b8SKris Buschelman 
131*1677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
132*1677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
133*1677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
134*1677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
135*1677d0b8SKris Buschelman 
136*1677d0b8SKris Buschelman     lu->flg = SAME_NONZERO_PATTERN;
137*1677d0b8SKris Buschelman   }
138*1677d0b8SKris Buschelman 
139*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
140*1677d0b8SKris Buschelman }
141*1677d0b8SKris Buschelman 
142*1677d0b8SKris Buschelman /*
143*1677d0b8SKris Buschelman    Note the r permutation is ignored
144*1677d0b8SKris Buschelman */
145*1677d0b8SKris Buschelman #undef __FUNCT__
146*1677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_UMFPACK"
147*1677d0b8SKris Buschelman int MatLUFactorSymbolic_SeqAIJ_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
148*1677d0b8SKris Buschelman {
149*1677d0b8SKris Buschelman   Mat_SeqAIJ          *mat=(Mat_SeqAIJ*)A->data;
150*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK  *lu;
151*1677d0b8SKris Buschelman   int                 ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ca;
152*1677d0b8SKris Buschelman   PetscScalar         *av=mat->a;
153*1677d0b8SKris Buschelman 
154*1677d0b8SKris Buschelman   PetscFunctionBegin;
155*1677d0b8SKris Buschelman   /* Create the factorization matrix F */
156*1677d0b8SKris Buschelman   ierr = MatCreateSeqAIJ(A->comm,m,n,PETSC_NULL,PETSC_NULL,F);CHKERRQ(ierr);
157*1677d0b8SKris Buschelman 
158*1677d0b8SKris Buschelman   (*F)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_UMFPACK;
159*1677d0b8SKris Buschelman   (*F)->ops->solve           = MatSolve_SeqAIJ_UMFPACK;
160*1677d0b8SKris Buschelman   (*F)->ops->destroy         = MatDestroy_SeqAIJ_UMFPACK;
161*1677d0b8SKris Buschelman   (*F)->factor               = FACTOR_LU;
162*1677d0b8SKris Buschelman   (*F)->assembled            = PETSC_TRUE;  /* required by -sles_view */
163*1677d0b8SKris Buschelman 
164*1677d0b8SKris Buschelman   ierr        = PetscNew(Mat_SeqAIJ_UMFPACK,&lu);CHKERRQ(ierr);
165*1677d0b8SKris Buschelman   (*F)->spptr = (void*)lu;
166*1677d0b8SKris Buschelman 
167*1677d0b8SKris Buschelman   /* initializations */
168*1677d0b8SKris Buschelman   /* ------------------------------------------------*/
169*1677d0b8SKris Buschelman   /* get the default control parameters */
170*1677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control) ;
171*1677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
172*1677d0b8SKris Buschelman 
173*1677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
174*1677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
175*1677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
176*1677d0b8SKris Buschelman 
177*1677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
178*1677d0b8SKris Buschelman   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);
179*1677d0b8SKris Buschelman   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);
180*1677d0b8SKris Buschelman   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);
181*1677d0b8SKris Buschelman 
182*1677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
183*1677d0b8SKris Buschelman   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);
184*1677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_relaxed_amalgamation","Control[UMFPACK_RELAXED_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED_AMALGAMATION],&lu->Control[UMFPACK_RELAXED_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr);
185*1677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_relaxed2_amalgamation","Control[UMFPACK_RELAXED2_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED2_AMALGAMATION],&lu->Control[UMFPACK_RELAXED2_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr);
186*1677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_relaxed3_amalgamation","Control[UMFPACK_RELAXED3_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED3_AMALGAMATION],&lu->Control[UMFPACK_RELAXED3_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr);
187*1677d0b8SKris Buschelman   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);
188*1677d0b8SKris Buschelman 
189*1677d0b8SKris Buschelman   /* Control parameters used by solve */
190*1677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
191*1677d0b8SKris Buschelman 
192*1677d0b8SKris Buschelman   /* use Petsc mat ordering (notice size is for the transpose) */
193*1677d0b8SKris Buschelman   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
194*1677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
195*1677d0b8SKris Buschelman     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
196*1677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
197*1677d0b8SKris Buschelman     ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr);
198*1677d0b8SKris Buschelman     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
199*1677d0b8SKris Buschelman   }
200*1677d0b8SKris Buschelman   PetscOptionsEnd();
201*1677d0b8SKris Buschelman 
202*1677d0b8SKris Buschelman   /* print the control parameters */
203*1677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
204*1677d0b8SKris Buschelman 
205*1677d0b8SKris Buschelman   /* symbolic factorization of A' */
206*1677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
207*1677d0b8SKris Buschelman   status = umfpack_di_qsymbolic(n,m,ai,aj,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ;
208*1677d0b8SKris Buschelman   if (status < 0){
209*1677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info) ;
210*1677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
211*1677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_symbolic failed");
212*1677d0b8SKris Buschelman   }
213*1677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
214*1677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ;
215*1677d0b8SKris Buschelman 
216*1677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
217*1677d0b8SKris Buschelman   lu->ai  = ai;
218*1677d0b8SKris Buschelman   lu->aj  = aj;
219*1677d0b8SKris Buschelman   lu->av  = av;
220*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
221*1677d0b8SKris Buschelman }
222*1677d0b8SKris Buschelman 
223*1677d0b8SKris Buschelman #undef __FUNCT__
224*1677d0b8SKris Buschelman #define __FUNCT__ "MatUseUMFPACK_SeqAIJ"
225*1677d0b8SKris Buschelman int MatUseUMFPACK_SeqAIJ(Mat A)
226*1677d0b8SKris Buschelman {
227*1677d0b8SKris Buschelman   PetscFunctionBegin;
228*1677d0b8SKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_UMFPACK;
229*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
230*1677d0b8SKris Buschelman }
231*1677d0b8SKris Buschelman 
232*1677d0b8SKris Buschelman /* used by -sles_view */
233*1677d0b8SKris Buschelman #undef __FUNCT__
234*1677d0b8SKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_UMFPACK"
235*1677d0b8SKris Buschelman int MatSeqAIJFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
236*1677d0b8SKris Buschelman {
237*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK      *lu= (Mat_SeqAIJ_UMFPACK*)A->spptr;
238*1677d0b8SKris Buschelman   int                     ierr;
239*1677d0b8SKris Buschelman   PetscFunctionBegin;
240*1677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
241*1677d0b8SKris Buschelman   if (A->ops->solve != MatSolve_SeqAIJ_UMFPACK) PetscFunctionReturn(0);
242*1677d0b8SKris Buschelman 
243*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
244*1677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
245*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
246*1677d0b8SKris Buschelman 
247*1677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
248*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
249*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
250*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
251*1677d0b8SKris Buschelman 
252*1677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
253*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
254*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED_AMALGAMATION]);CHKERRQ(ierr);
255*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED2_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED2_AMALGAMATION]);CHKERRQ(ierr);
256*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED3_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED3_AMALGAMATION]);CHKERRQ(ierr);
257*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
258*1677d0b8SKris Buschelman 
259*1677d0b8SKris Buschelman   /* Control parameters used by solve */
260*1677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
261*1677d0b8SKris Buschelman 
262*1677d0b8SKris Buschelman   /* mat ordering */
263*1677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
264*1677d0b8SKris Buschelman 
265*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
266*1677d0b8SKris Buschelman }
267*1677d0b8SKris Buschelman 
268*1677d0b8SKris Buschelman EXTERN_C_BEGIN
269*1677d0b8SKris Buschelman #undef __FUNCT__
270*1677d0b8SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_UMFPACK"
271*1677d0b8SKris Buschelman int MatCreate_SeqAIJ_UMFPACK(Mat A) {
272*1677d0b8SKris Buschelman   int                ierr;
273*1677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu;
274*1677d0b8SKris Buschelman 
275*1677d0b8SKris Buschelman   PetscFunctionBegin;
276*1677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
277*1677d0b8SKris Buschelman   ierr = MatUseUMFPACK_SeqAIJ(A);CHKERRQ(ierr);
278*1677d0b8SKris Buschelman 
279*1677d0b8SKris Buschelman   ierr                = PetscNew(Mat_SeqAIJ_UMFPACK,&lu);CHKERRQ(ierr);
280*1677d0b8SKris Buschelman   lu->MatView         = A->ops->view;
281*1677d0b8SKris Buschelman   lu->MatAssemblyEnd  = A->ops->assemblyend;
282*1677d0b8SKris Buschelman   lu->MatDestroy      = A->ops->destroy;
283*1677d0b8SKris Buschelman   lu->CleanUpUMFPACK  = PETSC_FALSE;
284*1677d0b8SKris Buschelman   A->spptr            = (void*)lu;
285*1677d0b8SKris Buschelman   A->ops->view        = MatView_SeqAIJ_UMFPACK;
286*1677d0b8SKris Buschelman   A->ops->assemblyend = MatAssemblyEnd_SeqAIJ_UMFPACK;
287*1677d0b8SKris Buschelman   A->ops->destroy     = MatDestroy_SeqAIJ_UMFPACK;
288*1677d0b8SKris Buschelman   PetscFunctionReturn(0);
289*1677d0b8SKris Buschelman }
290*1677d0b8SKris Buschelman EXTERN_C_END
291