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