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