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