1 2 /* 3 Provides an interface to the UMFPACKv5.1 sparse solver 4 5 When build with PETSC_USE_64BIT_INDICES this will use UF_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 UMFPACK UL_Long version MUST be built with 64 bit integers when used. 9 10 */ 11 #include <../src/mat/impls/aij/seq/aij.h> 12 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 #define umfpack_UMF_wsolve umfpack_zl_wsolve 18 #define umfpack_UMF_numeric umfpack_zl_numeric 19 #define umfpack_UMF_report_numeric umfpack_zl_report_numeric 20 #define umfpack_UMF_report_control umfpack_zl_report_control 21 #define umfpack_UMF_report_status umfpack_zl_report_status 22 #define umfpack_UMF_report_info umfpack_zl_report_info 23 #define umfpack_UMF_report_symbolic umfpack_zl_report_symbolic 24 #define umfpack_UMF_qsymbolic umfpack_zl_qsymbolic 25 #define umfpack_UMF_symbolic umfpack_zl_symbolic 26 #define umfpack_UMF_defaults umfpack_zl_defaults 27 28 #else 29 #define umfpack_UMF_free_symbolic umfpack_dl_free_symbolic 30 #define umfpack_UMF_free_numeric umfpack_dl_free_numeric 31 #define umfpack_UMF_wsolve umfpack_dl_wsolve 32 #define umfpack_UMF_numeric umfpack_dl_numeric 33 #define umfpack_UMF_report_numeric umfpack_dl_report_numeric 34 #define umfpack_UMF_report_control umfpack_dl_report_control 35 #define umfpack_UMF_report_status umfpack_dl_report_status 36 #define umfpack_UMF_report_info umfpack_dl_report_info 37 #define umfpack_UMF_report_symbolic umfpack_dl_report_symbolic 38 #define umfpack_UMF_qsymbolic umfpack_dl_qsymbolic 39 #define umfpack_UMF_symbolic umfpack_dl_symbolic 40 #define umfpack_UMF_defaults umfpack_dl_defaults 41 #endif 42 43 #else 44 #if defined(PETSC_USE_COMPLEX) 45 #define umfpack_UMF_free_symbolic umfpack_zi_free_symbolic 46 #define umfpack_UMF_free_numeric umfpack_zi_free_numeric 47 #define umfpack_UMF_wsolve umfpack_zi_wsolve 48 #define umfpack_UMF_numeric umfpack_zi_numeric 49 #define umfpack_UMF_report_numeric umfpack_zi_report_numeric 50 #define umfpack_UMF_report_control umfpack_zi_report_control 51 #define umfpack_UMF_report_status umfpack_zi_report_status 52 #define umfpack_UMF_report_info umfpack_zi_report_info 53 #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic 54 #define umfpack_UMF_qsymbolic umfpack_zi_qsymbolic 55 #define umfpack_UMF_symbolic umfpack_zi_symbolic 56 #define umfpack_UMF_defaults umfpack_zi_defaults 57 58 #else 59 #define umfpack_UMF_free_symbolic umfpack_di_free_symbolic 60 #define umfpack_UMF_free_numeric umfpack_di_free_numeric 61 #define umfpack_UMF_wsolve umfpack_di_wsolve 62 #define umfpack_UMF_numeric umfpack_di_numeric 63 #define umfpack_UMF_report_numeric umfpack_di_report_numeric 64 #define umfpack_UMF_report_control umfpack_di_report_control 65 #define umfpack_UMF_report_status umfpack_di_report_status 66 #define umfpack_UMF_report_info umfpack_di_report_info 67 #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic 68 #define umfpack_UMF_qsymbolic umfpack_di_qsymbolic 69 #define umfpack_UMF_symbolic umfpack_di_symbolic 70 #define umfpack_UMF_defaults umfpack_di_defaults 71 #endif 72 #endif 73 74 75 #define UF_long long long 76 #define UF_long_max LONG_LONG_MAX 77 #define UF_long_id "%lld" 78 79 EXTERN_C_BEGIN 80 #include <umfpack.h> 81 EXTERN_C_END 82 83 static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0}; 84 85 typedef struct { 86 void *Symbolic, *Numeric; 87 double Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W; 88 PetscInt *Wi,*ai,*aj,*perm_c; 89 PetscScalar *av; 90 MatStructure flg; 91 PetscBool PetscMatOrdering; 92 93 /* Flag to clean up UMFPACK objects during Destroy */ 94 PetscBool CleanUpUMFPACK; 95 } Mat_UMFPACK; 96 97 #undef __FUNCT__ 98 #define __FUNCT__ "MatDestroy_UMFPACK" 99 static PetscErrorCode MatDestroy_UMFPACK(Mat A) 100 { 101 PetscErrorCode ierr; 102 Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 103 104 PetscFunctionBegin; 105 if (lu->CleanUpUMFPACK) { 106 umfpack_UMF_free_symbolic(&lu->Symbolic); 107 umfpack_UMF_free_numeric(&lu->Numeric); 108 ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 109 ierr = PetscFree(lu->W);CHKERRQ(ierr); 110 if (lu->PetscMatOrdering) { 111 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 112 } 113 } 114 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatSolve_UMFPACK_Private" 120 static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag) 121 { 122 Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr; 123 PetscScalar *av=lu->av,*ba,*xa; 124 PetscErrorCode ierr; 125 PetscInt *ai=lu->ai,*aj=lu->aj,status; 126 127 PetscFunctionBegin; 128 /* solve Ax = b by umfpack_*_wsolve */ 129 /* ----------------------------------*/ 130 131 ierr = VecGetArray(b,&ba); 132 ierr = VecGetArray(x,&xa); 133 #if defined(PETSC_USE_COMPLEX) 134 status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL, 135 lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 136 #else 137 status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 138 #endif 139 umfpack_UMF_report_info(lu->Control, lu->Info); 140 if (status < 0){ 141 umfpack_UMF_report_status(lu->Control, status); 142 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed"); 143 } 144 145 ierr = VecRestoreArray(b,&ba);CHKERRQ(ierr); 146 ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149 150 #undef __FUNCT__ 151 #define __FUNCT__ "MatSolve_UMFPACK" 152 static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x) 153 { 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 158 ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 #undef __FUNCT__ 163 #define __FUNCT__ "MatSolveTranspose_UMFPACK" 164 static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x) 165 { 166 PetscErrorCode ierr; 167 168 PetscFunctionBegin; 169 /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 170 ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr); 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK" 176 static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info) 177 { 178 Mat_UMFPACK *lu=(Mat_UMFPACK*)(F)->spptr; 179 PetscErrorCode ierr; 180 PetscInt *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status; 181 PetscScalar *av=lu->av; 182 183 PetscFunctionBegin; 184 /* numeric factorization of A' */ 185 /* ----------------------------*/ 186 187 if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){ 188 umfpack_UMF_free_numeric(&lu->Numeric); 189 } 190 #if defined(PETSC_USE_COMPLEX) 191 status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 192 #else 193 status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 194 #endif 195 if (status < 0) { 196 umfpack_UMF_report_status(lu->Control, status); 197 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed"); 198 } 199 /* report numeric factorization of A' when Control[PRL] > 3 */ 200 (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control); 201 202 if (lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */ 203 /* allocate working space to be used by Solve */ 204 ierr = PetscMalloc(m * sizeof(PetscInt), &lu->Wi);CHKERRQ(ierr); 205 ierr = PetscMalloc(5*m * sizeof(PetscScalar), &lu->W);CHKERRQ(ierr); 206 } 207 208 lu->flg = SAME_NONZERO_PATTERN; 209 lu->CleanUpUMFPACK = PETSC_TRUE; 210 F->ops->solve = MatSolve_UMFPACK; 211 F->ops->solvetranspose = MatSolveTranspose_UMFPACK; 212 PetscFunctionReturn(0); 213 } 214 215 /* 216 Note the r permutation is ignored 217 */ 218 #undef __FUNCT__ 219 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK" 220 static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 221 { 222 Mat_SeqAIJ *mat=(Mat_SeqAIJ*)A->data; 223 Mat_UMFPACK *lu = (Mat_UMFPACK*)(F->spptr); 224 PetscErrorCode ierr; 225 PetscInt i,m=A->rmap->n,n=A->cmap->n; 226 const PetscInt *ra; 227 PetscInt status; 228 PetscScalar *av=mat->a; 229 230 PetscFunctionBegin; 231 if (lu->PetscMatOrdering) { 232 ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 233 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 234 /* we cannot simply memcpy on 64 bit archs */ 235 for(i = 0; i < m; i++) lu->perm_c[i] = ra[i]; 236 ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 237 } 238 239 lu->ai = mat->i; 240 lu->aj = mat->j; 241 242 /* print the control parameters */ 243 if(lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control); 244 245 /* symbolic factorization of A' */ 246 /* ---------------------------------------------------------------------- */ 247 if (lu->PetscMatOrdering) { /* use Petsc row ordering */ 248 #if !defined(PETSC_USE_COMPLEX) 249 status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 250 #else 251 status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,NULL,NULL, 252 lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 253 #endif 254 } else { /* use Umfpack col ordering */ 255 #if !defined(PETSC_USE_COMPLEX) 256 status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,av,&lu->Symbolic,lu->Control,lu->Info); 257 #else 258 status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info); 259 #endif 260 } 261 if (status < 0){ 262 umfpack_UMF_report_info(lu->Control, lu->Info); 263 umfpack_UMF_report_status(lu->Control, status); 264 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed"); 265 } 266 /* report sumbolic factorization of A' when Control[PRL] > 3 */ 267 (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control); 268 269 lu->flg = DIFFERENT_NONZERO_PATTERN; 270 lu->av = av; 271 lu->CleanUpUMFPACK = PETSC_TRUE; 272 (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatFactorInfo_UMFPACK" 278 static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) 279 { 280 Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 /* check if matrix is UMFPACK type */ 285 if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 286 287 ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 288 /* Control parameters used by reporting routiones */ 289 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 290 291 /* Control parameters used by symbolic factorization */ 292 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr); 293 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 294 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 295 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr); 296 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);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->PetscMatOrdering) { 312 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr); 313 } 314 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "MatView_UMFPACK" 320 static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer) 321 { 322 PetscErrorCode ierr; 323 PetscBool iascii; 324 PetscViewerFormat format; 325 326 PetscFunctionBegin; 327 ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr); 328 329 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 330 if (iascii) { 331 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 332 if (format == PETSC_VIEWER_ASCII_INFO) { 333 ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 334 } 335 } 336 PetscFunctionReturn(0); 337 } 338 339 EXTERN_C_BEGIN 340 #undef __FUNCT__ 341 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack" 342 PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type) 343 { 344 PetscFunctionBegin; 345 *type = MATSOLVERUMFPACK; 346 PetscFunctionReturn(0); 347 } 348 EXTERN_C_END 349 350 351 /*MC 352 MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 353 via the external package UMFPACK. 354 355 ./configure --download-umfpack to install PETSc to use UMFPACK 356 357 Consult UMFPACK documentation for more information about the Control parameters 358 which correspond to the options database keys below. 359 360 Options Database Keys: 361 + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 362 . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2 363 . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 364 . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW] 365 . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE] 366 . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 367 . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE] 368 . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ] 369 . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE] 370 . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 371 . -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE] 372 . -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX 373 . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 374 . -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL] 375 - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 376 377 Level: beginner 378 379 .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage 380 M*/ 381 EXTERN_C_BEGIN 382 #undef __FUNCT__ 383 #define __FUNCT__ "MatGetFactor_seqaij_umfpack" 384 PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F) 385 { 386 Mat B; 387 Mat_UMFPACK *lu; 388 PetscErrorCode ierr; 389 PetscInt m=A->rmap->n,n=A->cmap->n,idx; 390 391 const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"}, 392 *scale[]={"NONE","SUM","MAX"}; 393 PetscBool flg; 394 395 PetscFunctionBegin; 396 /* Create the factorization matrix F */ 397 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 398 ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 399 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 400 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 401 ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr); 402 B->spptr = lu; 403 B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 404 B->ops->destroy = MatDestroy_UMFPACK; 405 B->ops->view = MatView_UMFPACK; 406 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_umfpack",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr); 407 B->factortype = MAT_FACTOR_LU; 408 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 409 B->preallocated = PETSC_TRUE; 410 411 /* initializations */ 412 /* ------------------------------------------------*/ 413 /* get the default control parameters */ 414 umfpack_UMF_defaults(lu->Control); 415 lu->perm_c = PETSC_NULL; /* use defaul UMFPACK col permutation */ 416 lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */ 417 418 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 419 /* Control parameters used by reporting routiones */ 420 ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr); 421 422 /* Control parameters for symbolic factorization */ 423 ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr); 424 if (flg) { 425 switch (idx){ 426 case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break; 427 case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break; 428 case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break; 429 } 430 } 431 ierr = PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,sizeof UmfpackOrderingTypes/sizeof UmfpackOrderingTypes[0],UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg);CHKERRQ(ierr); 432 if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx; 433 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); 434 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); 435 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); 436 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); 437 ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr); 438 ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr); 439 440 /* Control parameters used by numeric factorization */ 441 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); 442 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); 443 ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 444 if (flg) { 445 switch (idx){ 446 case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break; 447 case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break; 448 case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break; 449 } 450 } 451 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); 452 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); 453 ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr); 454 455 /* Control parameters used by solve */ 456 ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr); 457 458 /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */ 459 ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOrdering);CHKERRQ(ierr); 460 PetscOptionsEnd(); 461 *F = B; 462 PetscFunctionReturn(0); 463 } 464 EXTERN_C_END 465 466