1 #define PETSCMAT_DLL 2 3 /* 4 Provides an interface to the UMFPACKv5.1 sparse solver 5 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. 10 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). 16 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. 20 21 We thank Christophe Geuzaine <geuzaine@acm.caltech.edu> for upgrading this interface to the UMFPACKv5.1 22 */ 23 24 #include "src/mat/impls/aij/seq/aij.h" 25 26 EXTERN_C_BEGIN 27 #include "umfpack.h" 28 EXTERN_C_END 29 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; 37 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; 48 49 EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*); 50 51 EXTERN_C_BEGIN 52 #undef __FUNCT__ 53 #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ" 54 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat) 55 { 56 PetscErrorCode ierr; 57 Mat B=*newmat; 58 Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 59 60 PetscFunctionBegin; 61 if (reuse == MAT_INITIAL_MATRIX) { 62 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 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; 70 71 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_umfpack_C","",PETSC_NULL);CHKERRQ(ierr); 72 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_umfpack_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 73 74 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 75 *newmat = B; 76 PetscFunctionReturn(0); 77 } 78 EXTERN_C_END 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "MatDestroy_UMFPACK" 82 PetscErrorCode MatDestroy_UMFPACK(Mat A) 83 { 84 PetscErrorCode ierr; 85 Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 86 87 PetscFunctionBegin; 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 ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 97 ierr = PetscFree(lu->W);CHKERRQ(ierr); 98 if(sizeof(UF_long) != sizeof(int)){ 99 ierr = PetscFree(lu->ai);CHKERRQ(ierr); 100 ierr = PetscFree(lu->aj);CHKERRQ(ierr); 101 } 102 if (lu->PetscMatOdering) { 103 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 104 } 105 } 106 ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 107 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 111 #undef __FUNCT__ 112 #define __FUNCT__ "MatSolve_UMFPACK" 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; 117 PetscErrorCode ierr; 118 UF_long *ai=lu->ai,*aj=lu->aj,status; 119 120 PetscFunctionBegin; 121 /* solve Ax = b by umfpack_*_wsolve */ 122 /* ----------------------------------*/ 123 124 #if defined(PETSC_USE_COMPLEX) 125 ierr = VecConjugate(b); 126 #endif 127 128 ierr = VecGetArray(b,&ba); 129 ierr = VecGetArray(x,&xa); 130 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 148 149 ierr = VecRestoreArray(b,&ba); 150 ierr = VecRestoreArray(x,&xa); 151 152 #if defined(PETSC_USE_COMPLEX) 153 ierr = VecConjugate(b); 154 ierr = VecConjugate(x); 155 #endif 156 157 PetscFunctionReturn(0); 158 } 159 160 #undef __FUNCT__ 161 #define __FUNCT__ "MatLUFactorNumeric_UMFPACK" 162 PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F) 163 { 164 Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr; 165 PetscErrorCode ierr; 166 UF_long *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status; 167 PetscScalar *av=lu->av; 168 169 PetscFunctionBegin; 170 /* numeric factorization of A' */ 171 /* ----------------------------*/ 172 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 190 191 if (lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */ 192 /* allocate working space to be used by Solve */ 193 ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr); 194 #if defined(PETSC_USE_COMPLEX) 195 ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr); 196 #else 197 ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr); 198 #endif 199 } 200 201 lu->flg = SAME_NONZERO_PATTERN; 202 lu->CleanUpUMFPACK = PETSC_TRUE; 203 PetscFunctionReturn(0); 204 } 205 206 /* 207 Note the r permutation is ignored 208 */ 209 #undef __FUNCT__ 210 #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK" 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; 216 PetscErrorCode ierr; 217 int i,m=A->rmap.n,n=A->cmap.n,*ra,idx; 218 UF_long status; 219 220 PetscScalar *av=mat->a; 221 const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"}, 222 *scale[]={"NONE","SUM","MAX"}; 223 PetscTruth flg; 224 225 PetscFunctionBegin; 226 /* Create the factorization matrix F */ 227 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 228 ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 229 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 230 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 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 */ 236 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 */ 249 250 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 251 /* Control parameters used by reporting routiones */ 252 ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr); 253 254 /* Control parameters for symbolic factorization */ 255 ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr); 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 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); 265 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); 266 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); 267 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); 268 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); 269 ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr); 270 ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr); 271 272 /* Control parameters used by numeric factorization */ 273 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); 274 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); 275 ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 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 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); 284 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); 285 ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr); 286 287 /* Control parameters used by solve */ 288 ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr); 289 290 /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */ 291 ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr); 292 if (lu->PetscMatOdering) { 293 ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 294 ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr); 295 /* we cannot simply memcpy on 64 bit archs */ 296 for(i = 0; i < m; i++) lu->perm_c[i] = ra[i]; 297 ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 298 } 299 PetscOptionsEnd(); 300 301 if(sizeof(UF_long) != sizeof(int)){ 302 /* we cannot directly use mat->i and mat->j on 64 bit archs */ 303 ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr); 304 ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr); 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 } 312 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 319 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 353 354 lu->flg = DIFFERENT_NONZERO_PATTERN; 355 lu->av = av; 356 lu->CleanUpUMFPACK = PETSC_TRUE; 357 *F = B; 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "MatAssemblyEnd_UMFPACK" 363 PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) 364 { 365 PetscErrorCode ierr; 366 Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 367 368 PetscFunctionBegin; 369 ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 370 lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 371 A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 372 PetscFunctionReturn(0); 373 } 374 375 /* used by -ksp_view */ 376 #undef __FUNCT__ 377 #define __FUNCT__ "MatFactorInfo_UMFPACK" 378 PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) 379 { 380 Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr; 381 PetscErrorCode ierr; 382 383 PetscFunctionBegin; 384 /* check if matrix is UMFPACK type */ 385 if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 386 387 ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 388 /* Control parameters used by reporting routiones */ 389 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 390 391 /* Control parameters used by symbolic factorization */ 392 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr); 393 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 394 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 395 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr); 396 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 397 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr); 398 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr); 399 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr); 400 401 /* Control parameters used by numeric factorization */ 402 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr); 403 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr); 404 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr); 405 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr); 406 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr); 407 408 /* Control parameters used by solve */ 409 ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr); 410 411 /* mat ordering */ 412 if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer," UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr); 413 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "MatView_UMFPACK" 419 PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer) 420 { 421 PetscErrorCode ierr; 422 PetscTruth iascii; 423 PetscViewerFormat format; 424 Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 425 426 PetscFunctionBegin; 427 ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 428 429 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 430 if (iascii) { 431 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 432 if (format == PETSC_VIEWER_ASCII_INFO) { 433 ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 434 } 435 } 436 PetscFunctionReturn(0); 437 } 438 439 EXTERN_C_BEGIN 440 #undef __FUNCT__ 441 #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK" 442 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat) 443 { 444 PetscErrorCode ierr; 445 Mat B=*newmat; 446 Mat_UMFPACK *lu; 447 448 PetscFunctionBegin; 449 if (reuse == MAT_INITIAL_MATRIX) { 450 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 451 } 452 453 ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr); 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; 460 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; 467 468 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C", 469 "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr); 470 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C", 471 "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr); 472 ierr = PetscInfo(A,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr); 473 ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr); 474 *newmat = B; 475 PetscFunctionReturn(0); 476 } 477 EXTERN_C_END 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "MatDuplicate_UMFPACK" 481 PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) 482 { 483 PetscErrorCode ierr; 484 Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 485 486 PetscFunctionBegin; 487 ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 488 ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 /*MC 493 MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 494 via the external package UMFPACK. 495 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). 500 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. 504 505 Consult UMFPACK documentation for more information about the Control parameters 506 which correspond to the options database keys below. 507 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] 516 517 Level: beginner 518 519 .seealso: PCLU 520 M*/ 521 522 EXTERN_C_BEGIN 523 #undef __FUNCT__ 524 #define __FUNCT__ "MatCreate_UMFPACK" 525 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A) 526 { 527 PetscErrorCode ierr; 528 529 PetscFunctionBegin; 530 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 531 ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 EXTERN_C_END 535 536 537 538