Actual source code: umfpack.c


  2: /*
  3:    Provides an interface to the UMFPACK sparse solver available through SuiteSparse version 4.2.1

  5:    When build with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
  6:    integer type in UMFPACK, otherwise it will use int. This means
  7:    all integers in this file as simply declared as PetscInt. Also it means
  8:    that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit only]

 10: */
 11: #include <../src/mat/impls/aij/seq/aij.h>

 13: #if defined(PETSC_USE_64BIT_INDICES)
 14: #if defined(PETSC_USE_COMPLEX)
 15: #define umfpack_UMF_free_symbolic                      umfpack_zl_free_symbolic
 16: #define umfpack_UMF_free_numeric                       umfpack_zl_free_numeric
 17: /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
 18: #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k,l,m,n) umfpack_zl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k,l,(SuiteSparse_long*)m,n)
 19: #define umfpack_UMF_numeric(a,b,c,d,e,f,g,h)          umfpack_zl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g,h)
 20: #define umfpack_UMF_report_numeric                    umfpack_zl_report_numeric
 21: #define umfpack_UMF_report_control                    umfpack_zl_report_control
 22: #define umfpack_UMF_report_status                     umfpack_zl_report_status
 23: #define umfpack_UMF_report_info                       umfpack_zl_report_info
 24: #define umfpack_UMF_report_symbolic                   umfpack_zl_report_symbolic
 25: #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i,j)    umfpack_zl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,(SuiteSparse_long*)g,h,i,j)
 26: #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h,i)       umfpack_zl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h,i)
 27: #define umfpack_UMF_defaults                          umfpack_zl_defaults

 29: #else
 30: #define umfpack_UMF_free_symbolic                  umfpack_dl_free_symbolic
 31: #define umfpack_UMF_free_numeric                   umfpack_dl_free_numeric
 32: #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k)  umfpack_dl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k)
 33: #define umfpack_UMF_numeric(a,b,c,d,e,f,g)         umfpack_dl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g)
 34: #define umfpack_UMF_report_numeric                 umfpack_dl_report_numeric
 35: #define umfpack_UMF_report_control                 umfpack_dl_report_control
 36: #define umfpack_UMF_report_status                  umfpack_dl_report_status
 37: #define umfpack_UMF_report_info                    umfpack_dl_report_info
 38: #define umfpack_UMF_report_symbolic                umfpack_dl_report_symbolic
 39: #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i)   umfpack_dl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,(SuiteSparse_long*)f,g,h,i)
 40: #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h)      umfpack_dl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h)
 41: #define umfpack_UMF_defaults                       umfpack_dl_defaults
 42: #endif

 44: #else
 45: #if defined(PETSC_USE_COMPLEX)
 46: #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
 47: #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
 48: #define umfpack_UMF_wsolve          umfpack_zi_wsolve
 49: #define umfpack_UMF_numeric         umfpack_zi_numeric
 50: #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
 51: #define umfpack_UMF_report_control  umfpack_zi_report_control
 52: #define umfpack_UMF_report_status   umfpack_zi_report_status
 53: #define umfpack_UMF_report_info     umfpack_zi_report_info
 54: #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
 55: #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
 56: #define umfpack_UMF_symbolic        umfpack_zi_symbolic
 57: #define umfpack_UMF_defaults        umfpack_zi_defaults

 59: #else
 60: #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
 61: #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
 62: #define umfpack_UMF_wsolve          umfpack_di_wsolve
 63: #define umfpack_UMF_numeric         umfpack_di_numeric
 64: #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
 65: #define umfpack_UMF_report_control  umfpack_di_report_control
 66: #define umfpack_UMF_report_status   umfpack_di_report_status
 67: #define umfpack_UMF_report_info     umfpack_di_report_info
 68: #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
 69: #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
 70: #define umfpack_UMF_symbolic        umfpack_di_symbolic
 71: #define umfpack_UMF_defaults        umfpack_di_defaults
 72: #endif
 73: #endif

 75: EXTERN_C_BEGIN
 76: #include <umfpack.h>
 77: EXTERN_C_END

 79: static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};

 81: typedef struct {
 82:   void         *Symbolic, *Numeric;
 83:   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
 84:   PetscInt     *Wi,*perm_c;
 85:   Mat          A;               /* Matrix used for factorization */
 86:   MatStructure flg;

 88:   /* Flag to clean up UMFPACK objects during Destroy */
 89:   PetscBool CleanUpUMFPACK;
 90: } Mat_UMFPACK;

 92: static PetscErrorCode MatDestroy_UMFPACK(Mat A)
 93: {
 94:   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->data;

 96:   if (lu->CleanUpUMFPACK) {
 97:     umfpack_UMF_free_symbolic(&lu->Symbolic);
 98:     umfpack_UMF_free_numeric(&lu->Numeric);
 99:     PetscFree(lu->Wi);
100:     PetscFree(lu->W);
101:     PetscFree(lu->perm_c);
102:   }
103:   MatDestroy(&lu->A);
104:   PetscFree(A->data);
105:   return 0;
106: }

108: static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
109: {
110:   Mat_UMFPACK       *lu = (Mat_UMFPACK*)A->data;
111:   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)lu->A->data;
112:   PetscScalar       *av = a->a,*xa;
113:   const PetscScalar *ba;
114:   PetscInt          *ai = a->i,*aj = a->j,status;
115:   static PetscBool  cite = PETSC_FALSE;

117:   if (!A->rmap->n) return 0;
118:   PetscCitationsRegister("@article{davis2004algorithm,\n  title={Algorithm 832: {UMFPACK} V4.3---An Unsymmetric-Pattern Multifrontal Method},\n  author={Davis, Timothy A},\n  journal={ACM Transactions on Mathematical Software (TOMS)},\n  volume={30},\n  number={2},\n  pages={196--199},\n  year={2004},\n  publisher={ACM}\n}\n",&cite);
119:   /* solve Ax = b by umfpack_*_wsolve */
120:   /* ----------------------------------*/

122:   if (!lu->Wi) {  /* first time, allocate working space for wsolve */
123:     PetscMalloc1(A->rmap->n,&lu->Wi);
124:     PetscMalloc1(5*A->rmap->n,&lu->W);
125:   }

127:   VecGetArrayRead(b,&ba);
128:   VecGetArray(x,&xa);
129: #if defined(PETSC_USE_COMPLEX)
130:   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
131: #else
132:   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
133: #endif
134:   umfpack_UMF_report_info(lu->Control, lu->Info);
135:   if (status < 0) {
136:     umfpack_UMF_report_status(lu->Control, status);
137:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
138:   }

140:   VecRestoreArrayRead(b,&ba);
141:   VecRestoreArray(x,&xa);
142:   return 0;
143: }

145: static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
146: {
147:   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
148:   MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);
149:   return 0;
150: }

152: static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
153: {
154:   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
155:   MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);
156:   return 0;
157: }

159: static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
160: {
161:   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F)->data;
162:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
163:   PetscInt       *ai = a->i,*aj=a->j,status;
164:   PetscScalar    *av = a->a;

166:   if (!A->rmap->n) return 0;
167:   /* numeric factorization of A' */
168:   /* ----------------------------*/

170:   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
171:     umfpack_UMF_free_numeric(&lu->Numeric);
172:   }
173: #if defined(PETSC_USE_COMPLEX)
174:   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
175: #else
176:   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
177: #endif
178:   if (status < 0) {
179:     umfpack_UMF_report_status(lu->Control, status);
180:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
181:   }
182:   /* report numeric factorization of A' when Control[PRL] > 3 */
183:   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);

185:   PetscObjectReference((PetscObject)A);
186:   MatDestroy(&lu->A);

188:   lu->A                  = A;
189:   lu->flg                = SAME_NONZERO_PATTERN;
190:   lu->CleanUpUMFPACK     = PETSC_TRUE;
191:   F->ops->solve          = MatSolve_UMFPACK;
192:   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
193:   return 0;
194: }

196: static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
197: {
198:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
199:   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->data);
200:   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
201: #if !defined(PETSC_USE_COMPLEX)
202:   PetscScalar    *av = a->a;
203: #endif
204:   const PetscInt *ra;
205:   PetscInt       status;

207:   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
208:   if (!n) return 0;
209:   if (r) {
210:     ISGetIndices(r,&ra);
211:     PetscMalloc1(m,&lu->perm_c);
212:     /* we cannot simply memcpy on 64 bit archs */
213:     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
214:     ISRestoreIndices(r,&ra);
215:   }

217:   /* print the control parameters */
218:   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);

220:   /* symbolic factorization of A' */
221:   /* ---------------------------------------------------------------------- */
222:   if (r) { /* use Petsc row ordering */
223: #if !defined(PETSC_USE_COMPLEX)
224:     status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
225: #else
226:     status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
227: #endif
228:   } else { /* use Umfpack col ordering */
229: #if !defined(PETSC_USE_COMPLEX)
230:     status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
231: #else
232:     status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
233: #endif
234:   }
235:   if (status < 0) {
236:     umfpack_UMF_report_info(lu->Control, lu->Info);
237:     umfpack_UMF_report_status(lu->Control, status);
238:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
239:   }
240:   /* report sumbolic factorization of A' when Control[PRL] > 3 */
241:   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);

243:   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
244:   lu->CleanUpUMFPACK        = PETSC_TRUE;
245:   return 0;
246: }

248: static PetscErrorCode MatView_Info_UMFPACK(Mat A,PetscViewer viewer)
249: {
250:   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->data;

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

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

259:   /* Control parameters used by symbolic factorization */
260:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);
261:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);
262:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);
263:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);
264:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);
265:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);
266:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);

268:   /* Control parameters used by numeric factorization */
269:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);
270:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);
271:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);
272:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);
273:   PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);

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

278:   /* mat ordering */
279:   if (!lu->perm_c) {
280:     PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);
281:   }
282:   return 0;
283: }

285: static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
286: {
287:   PetscBool         iascii;
288:   PetscViewerFormat format;

290:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
291:   if (iascii) {
292:     PetscViewerGetFormat(viewer,&format);
293:     if (format == PETSC_VIEWER_ASCII_INFO) {
294:       MatView_Info_UMFPACK(A,viewer);
295:     }
296:   }
297:   return 0;
298: }

300: PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A,MatSolverType *type)
301: {
302:   *type = MATSOLVERUMFPACK;
303:   return 0;
304: }

306: /*MC
307:   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
308:   via the external package UMFPACK.

310:   Use ./configure --download-suitesparse to install PETSc to use UMFPACK

312:   Use -pc_type lu -pc_factor_mat_solver_type umfpack to use this direct solver

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

317:   Options Database Keys:
318: + -mat_umfpack_ordering                - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE
319: . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
320: . -mat_umfpack_strategy <AUTO>         - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
321: . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
322: . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
323: . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
324: . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
325: . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
326: . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
327: . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
328: . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
329: . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
330: . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
331: . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
332: . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
333: - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]

335:    Level: beginner

337:    Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html

339: .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverType(), MatSolverType
340: M*/

342: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
343: {
344:   Mat            B;
345:   Mat_UMFPACK    *lu;
347:   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
348:   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"};
349:   const char     *scale[]   ={"NONE","SUM","MAX"};
350:   PetscBool      flg;

352:   /* Create the factorization matrix F */
353:   MatCreate(PetscObjectComm((PetscObject)A),&B);
354:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
355:   PetscStrallocpy("umfpack",&((PetscObject)B)->type_name);
356:   MatSetUp(B);

358:   PetscNewLog(B,&lu);

360:   B->data                   = lu;
361:   B->ops->getinfo          = MatGetInfo_External;
362:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
363:   B->ops->destroy          = MatDestroy_UMFPACK;
364:   B->ops->view             = MatView_UMFPACK;
365:   B->ops->matsolve         = NULL;

367:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_umfpack);

369:   B->factortype   = MAT_FACTOR_LU;
370:   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
371:   B->preallocated = PETSC_TRUE;

373:   PetscFree(B->solvertype);
374:   PetscStrallocpy(MATSOLVERUMFPACK,&B->solvertype);
375:   B->canuseordering = PETSC_TRUE;
376:   PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]);

378:   /* initializations */
379:   /* ------------------------------------------------*/
380:   /* get the default control parameters */
381:   umfpack_UMF_defaults(lu->Control);
382:   lu->perm_c                  = NULL; /* use defaul UMFPACK col permutation */
383:   lu->Control[UMFPACK_IRSTEP] = 0;          /* max num of iterative refinement steps to attempt */

385:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");
386:   /* Control parameters used by reporting routiones */
387:   PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL);

389:   /* Control parameters for symbolic factorization */
390:   PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);
391:   if (flg) {
392:     switch (idx) {
393:     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
394:     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
395:     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
396:     }
397:   }
398:   PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,sizeof(UmfpackOrderingTypes)/sizeof(UmfpackOrderingTypes[0]),UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg);
399:   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
400:   PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],NULL);
401:   PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],NULL);
402:   PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],NULL);
403:   PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],NULL);
404:   PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL);
405:   PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL);

407:   /* Control parameters used by numeric factorization */
408:   PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],NULL);
409:   PetscOptionsReal("-mat_umfpack_sym_pivot_tolerance","Control[UMFPACK_SYM_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],&lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],NULL);
410:   PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);
411:   if (flg) {
412:     switch (idx) {
413:     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
414:     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
415:     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
416:     }
417:   }
418:   PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);
419:   PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);
420:   PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL);

422:   /* Control parameters used by solve */
423:   PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL);
424:   PetscOptionsEnd();
425:   *F = B;
426:   return 0;
427: }

429: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
430: PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
431: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
432: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat,MatFactorType,Mat*);

434: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void)
435: {
436:   MatSolverTypeRegister(MATSOLVERUMFPACK,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack);
437:   MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQAIJ,  MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod);
438:   MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod);
439:   MatSolverTypeRegister(MATSOLVERKLU,MATSEQAIJ,      MAT_FACTOR_LU,MatGetFactor_seqaij_klu);
440:   MatSolverTypeRegister(MATSOLVERSPQR,MATSEQAIJ,     MAT_FACTOR_QR,MatGetFactor_seqaij_spqr);
441:   if (!PetscDefined(USE_COMPLEX)) {
442:     MatSolverTypeRegister(MATSOLVERSPQR,MATNORMAL,   MAT_FACTOR_QR,MatGetFactor_seqaij_spqr);
443:   }
444:   MatSolverTypeRegister(MATSOLVERSPQR,MATNORMALHERMITIAN, MAT_FACTOR_QR,MatGetFactor_seqaij_spqr);
445:   return 0;
446: }