Actual source code: umfpack.c

  1: #define PETSCMAT_DLL

  3: /* 
  4:    Provides an interface to the UMFPACKv5.1 sparse solver

  6:    This interface uses the "UF_long version" of the UMFPACK API
  7:    (*_dl_* and *_zl_* routines, instead of *_di_* and *_zi_* routines)
  8:    so that UMFPACK can address more than 2Gb of memory on 64 bit
  9:    machines.

 11:    If sizeof(UF_long) == 32 the interface only allocates the memory
 12:    necessary for UMFPACK's working arrays (W, Wi and possibly
 13:    perm_c). If sizeof(UF_long) == 64, in addition to allocating the
 14:    working arrays, the interface also has to re-allocate the matrix
 15:    index arrays (ai and aj, which must be stored as UF_long).

 17:    The interface is implemented for both real and complex
 18:    arithmetic. Complex numbers are assumed to be in packed format,
 19:    which requires UMFPACK >= 4.4.

 21:    We thank Christophe Geuzaine <geuzaine@acm.caltech.edu> for upgrading this interface to the UMFPACKv5.1
 22: */

 24:  #include src/mat/impls/aij/seq/aij.h

 27: #include "umfpack.h"

 30: typedef struct {
 31:   void         *Symbolic, *Numeric;
 32:   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
 33:   UF_long      *Wi,*ai,*aj,*perm_c;
 34:   PetscScalar  *av;
 35:   MatStructure flg;
 36:   PetscTruth   PetscMatOdering;

 38:   /* A few function pointers for inheritance */
 39:   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
 40:   PetscErrorCode (*MatView)(Mat,PetscViewer);
 41:   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
 42:   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
 43:   PetscErrorCode (*MatDestroy)(Mat);
 44: 
 45:   /* Flag to clean up UMFPACK objects during Destroy */
 46:   PetscTruth CleanUpUMFPACK;
 47: } Mat_UMFPACK;

 49: EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);

 54: PetscErrorCode  MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
 55: {
 57:   Mat         B=*newmat;
 58:   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;

 61:   if (reuse == MAT_INITIAL_MATRIX) {
 62:     MatDuplicate(A,MAT_COPY_VALUES,&B);
 63:   }
 64:   /* Reset the original function pointers */
 65:   B->ops->duplicate        = lu->MatDuplicate;
 66:   B->ops->view             = lu->MatView;
 67:   B->ops->assemblyend      = lu->MatAssemblyEnd;
 68:   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
 69:   B->ops->destroy          = lu->MatDestroy;

 71:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_umfpack_C","",PETSC_NULL);
 72:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_umfpack_seqaij_C","",PETSC_NULL);

 74:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
 75:   *newmat = B;
 76:   return(0);
 77: }

 82: PetscErrorCode MatDestroy_UMFPACK(Mat A)
 83: {
 85:   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;

 88:   if (lu->CleanUpUMFPACK) {
 89: #if defined(PETSC_USE_COMPLEX)
 90:     umfpack_zl_free_symbolic(&lu->Symbolic);
 91:     umfpack_zl_free_numeric(&lu->Numeric);
 92: #else
 93:     umfpack_dl_free_symbolic(&lu->Symbolic);
 94:     umfpack_dl_free_numeric(&lu->Numeric);
 95: #endif
 96:     PetscFree(lu->Wi);
 97:     PetscFree(lu->W);
 98:     if(sizeof(UF_long) != sizeof(int)){
 99:       PetscFree(lu->ai);
100:       PetscFree(lu->aj);
101:     }
102:     if (lu->PetscMatOdering) {
103:       PetscFree(lu->perm_c);
104:     }
105:   }
106:   MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);
107:   (*A->ops->destroy)(A);
108:   return(0);
109: }

113: PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
114: {
115:   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
116:   PetscScalar *av=lu->av,*ba,*xa;
118:   UF_long     *ai=lu->ai,*aj=lu->aj,status;
119: 
121:   /* solve Ax = b by umfpack_*_wsolve */
122:   /* ----------------------------------*/

124: #if defined(PETSC_USE_COMPLEX)
125:   VecConjugate(b);
126: #endif

128:   VecGetArray(b,&ba);
129:   VecGetArray(x,&xa);

131: #if defined(PETSC_USE_COMPLEX)
132:   status = umfpack_zl_wsolve(UMFPACK_At,ai,aj,(double*)av,NULL,(double*)xa,NULL,(double*)ba,NULL,
133:                              lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
134:   umfpack_zl_report_info(lu->Control, lu->Info);
135:   if (status < 0){
136:     umfpack_zl_report_status(lu->Control, status);
137:     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_wsolve failed");
138:   }
139: #else
140:   status = umfpack_dl_wsolve(UMFPACK_At,ai,aj,av,xa,ba,
141:                              lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
142:   umfpack_dl_report_info(lu->Control, lu->Info);
143:   if (status < 0){
144:     umfpack_dl_report_status(lu->Control, status);
145:     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_wsolve failed");
146:   }
147: #endif

149:   VecRestoreArray(b,&ba);
150:   VecRestoreArray(x,&xa);

152: #if defined(PETSC_USE_COMPLEX)
153:   VecConjugate(b);
154:   VecConjugate(x);
155: #endif

157:   return(0);
158: }

162: PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
163: {
164:   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
166:   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status;
167:   PetscScalar *av=lu->av;

170:   /* numeric factorization of A' */
171:   /* ----------------------------*/

173: #if defined(PETSC_USE_COMPLEX)
174:   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
175:     umfpack_zl_free_numeric(&lu->Numeric);
176:   }
177:   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
178:   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
179:   /* report numeric factorization of A' when Control[PRL] > 3 */
180:   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
181: #else
182:   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
183:     umfpack_dl_free_numeric(&lu->Numeric);
184:   }
185:   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
186:   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
187:   /* report numeric factorization of A' when Control[PRL] > 3 */
188:   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
189: #endif

191:   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
192:     /* allocate working space to be used by Solve */
193:     PetscMalloc(m * sizeof(UF_long), &lu->Wi);
194: #if defined(PETSC_USE_COMPLEX)
195:     PetscMalloc(2*5*m * sizeof(double), &lu->W);
196: #else
197:     PetscMalloc(5*m * sizeof(double), &lu->W);
198: #endif
199:   }

201:   lu->flg = SAME_NONZERO_PATTERN;
202:   lu->CleanUpUMFPACK = PETSC_TRUE;
203:   return(0);
204: }

206: /*
207:    Note the r permutation is ignored
208: */
211: PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
212: {
213:   Mat         B;
214:   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
215:   Mat_UMFPACK *lu;
217:   int         i,m=A->rmap.n,n=A->cmap.n,*ra,idx;
218:   UF_long     status;

220:   PetscScalar *av=mat->a;
221:   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
222:               *scale[]={"NONE","SUM","MAX"};
223:   PetscTruth  flg;
224: 
226:   /* Create the factorization matrix F */
227:   MatCreate(((PetscObject)A)->comm,&B);
228:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
229:   MatSetType(B,((PetscObject)A)->type_name);
230:   MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
231: 
232:   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
233:   B->ops->solve           = MatSolve_UMFPACK;
234:   B->factor               = FACTOR_LU;
235:   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */

237:   lu = (Mat_UMFPACK*)(B->spptr);
238: 
239:   /* initializations */
240:   /* ------------------------------------------------*/
241:   /* get the default control parameters */
242: #if defined(PETSC_USE_COMPLEX)
243:   umfpack_zl_defaults(lu->Control);
244: #else
245:   umfpack_dl_defaults(lu->Control);
246: #endif
247:   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
248:   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */

250:   PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");
251:   /* Control parameters used by reporting routiones */
252:   PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);

254:   /* Control parameters for symbolic factorization */
255:   PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);
256:   if (flg) {
257:     switch (idx){
258:     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
259:     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
260:     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
261:     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
262:     }
263:   }
264:   PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],PETSC_NULL);
265:   PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],PETSC_NULL);
266:   PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],PETSC_NULL);
267:   PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],PETSC_NULL);
268:   PetscOptionsReal("-mat_umfpack_2by2_tolerance","Control[UMFPACK_2BY2_TOLERANCE]","None",lu->Control[UMFPACK_2BY2_TOLERANCE],&lu->Control[UMFPACK_2BY2_TOLERANCE],PETSC_NULL);
269:   PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);
270:   PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);

272:   /* Control parameters used by numeric factorization */
273:   PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],PETSC_NULL);
274:   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);
275:   PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);
276:   if (flg) {
277:     switch (idx){
278:     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
279:     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
280:     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
281:     }
282:   }
283:   PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);
284:   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);
285:   PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);

287:   /* Control parameters used by solve */
288:   PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);
289: 
290:   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
291:   PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);
292:   if (lu->PetscMatOdering) {
293:     ISGetIndices(r,&ra);
294:     PetscMalloc(m*sizeof(UF_long),&lu->perm_c);
295:     /* we cannot simply memcpy on 64 bit archs */
296:     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
297:     ISRestoreIndices(r,&ra);
298:   }
299:   PetscOptionsEnd();

301:   if(sizeof(UF_long) != sizeof(int)){
302:     /* we cannot directly use mat->i and mat->j on 64 bit archs */
303:     PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);
304:     PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);
305:     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
306:     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
307:   }
308:   else{
309:     lu->ai = (UF_long*)mat->i;
310:     lu->aj = (UF_long*)mat->j;
311:   }

313:   /* print the control parameters */
314: #if defined(PETSC_USE_COMPLEX)
315:   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
316: #else
317:   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
318: #endif

320:   /* symbolic factorization of A' */
321:   /* ---------------------------------------------------------------------- */
322: #if defined(PETSC_USE_COMPLEX)
323:   if (lu->PetscMatOdering) { /* use Petsc row ordering */
324:     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
325:                                   lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
326:   } else { /* use Umfpack col ordering */
327:     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
328:                                  &lu->Symbolic,lu->Control,lu->Info);
329:   }
330:   if (status < 0){
331:     umfpack_zl_report_info(lu->Control, lu->Info);
332:     umfpack_zl_report_status(lu->Control, status);
333:     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
334:   }
335:   /* report sumbolic factorization of A' when Control[PRL] > 3 */
336:   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
337: #else
338:   if (lu->PetscMatOdering) { /* use Petsc row ordering */
339:     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
340:                                   lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
341:   } else { /* use Umfpack col ordering */
342:     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
343:                                  &lu->Symbolic,lu->Control,lu->Info);
344:   }
345:   if (status < 0){
346:     umfpack_dl_report_info(lu->Control, lu->Info);
347:     umfpack_dl_report_status(lu->Control, status);
348:     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
349:   }
350:   /* report sumbolic factorization of A' when Control[PRL] > 3 */
351:   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
352: #endif

354:   lu->flg = DIFFERENT_NONZERO_PATTERN;
355:   lu->av  = av;
356:   lu->CleanUpUMFPACK = PETSC_TRUE;
357:   *F = B;
358:   return(0);
359: }

363: PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
364: {
366:   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);

369:   (*lu->MatAssemblyEnd)(A,mode);
370:   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
371:   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
372:   return(0);
373: }

375: /* used by -ksp_view */
378: PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
379: {
380:   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;

384:   /* check if matrix is UMFPACK type */
385:   if (A->ops->solve != MatSolve_UMFPACK) return(0);

387:   PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");
388:   /* Control parameters used by reporting routiones */
389:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);

391:   /* Control parameters used by symbolic factorization */
392:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);
393:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);
394:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);
395:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);
396:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);
397:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);
398:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);
399:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);

401:   /* Control parameters used by numeric factorization */
402:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);
403:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);
404:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);
405:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);
406:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);

408:   /* Control parameters used by solve */
409:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);

411:   /* mat ordering */
412:   if(!lu->PetscMatOdering) PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");

414:   return(0);
415: }

419: PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
420: {
422:   PetscTruth        iascii;
423:   PetscViewerFormat format;
424:   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);

427:   (*lu->MatView)(A,viewer);

429:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
430:   if (iascii) {
431:     PetscViewerGetFormat(viewer,&format);
432:     if (format == PETSC_VIEWER_ASCII_INFO) {
433:       MatFactorInfo_UMFPACK(A,viewer);
434:     }
435:   }
436:   return(0);
437: }

442: PetscErrorCode  MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
443: {
445:   Mat            B=*newmat;
446:   Mat_UMFPACK    *lu;

449:   if (reuse == MAT_INITIAL_MATRIX) {
450:     MatDuplicate(A,MAT_COPY_VALUES,&B);
451:   }

453:   PetscNewLog(B,Mat_UMFPACK,&lu);
454:   lu->MatDuplicate         = A->ops->duplicate;
455:   lu->MatView              = A->ops->view;
456:   lu->MatAssemblyEnd       = A->ops->assemblyend;
457:   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
458:   lu->MatDestroy           = A->ops->destroy;
459:   lu->CleanUpUMFPACK       = PETSC_FALSE;

461:   B->spptr                 = (void*)lu;
462:   B->ops->duplicate        = MatDuplicate_UMFPACK;
463:   B->ops->view             = MatView_UMFPACK;
464:   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
465:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
466:   B->ops->destroy          = MatDestroy_UMFPACK;

468:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
469:                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);
470:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
471:                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);
472:   PetscInfo(A,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");
473:   PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);
474:   *newmat = B;
475:   return(0);
476: }

481: PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
482: {
484:   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;

487:   (*lu->MatDuplicate)(A,op,M);
488:   PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));
489:   return(0);
490: }

492: /*MC
493:   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 
494:   via the external package UMFPACK.

496:   If UMFPACK is installed (see the manual for
497:   instructions on how to declare the existence of external packages),
498:   a matrix type can be constructed which invokes UMFPACK solvers.
499:   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).

501:   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is 
502:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
503:   the MATSEQAIJ type without data copy.

505:   Consult UMFPACK documentation for more information about the Control parameters
506:   which correspond to the options database keys below.

508:   Options Database Keys:
509: + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
510: . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
511: . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
512: . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
513: . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
514: . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
515: - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]

517:    Level: beginner

519: .seealso: PCLU
520: M*/

525: PetscErrorCode  MatCreate_UMFPACK(Mat A)
526: {

530:   MatSetType(A,MATSEQAIJ);
531:   MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);
532:   return(0);
533: }