1 2 /* -------------------------------------------------------------------- 3 4 This file implements a subclass of the SeqAIJ matrix class that uses 5 the SuperLU sparse solver. You can use this as a starting point for 6 implementing your own subclass of a PETSc matrix class. 7 8 This demonstrates a way to make an implementation inheritence of a PETSc 9 matrix type. This means constructing a new matrix type (SuperLU) by changing some 10 of the methods of the previous type (SeqAIJ), adding additional data, and possibly 11 additional method. (See any book on object oriented programming). 12 */ 13 14 /* 15 Defines the data structure for the base matrix type (SeqAIJ) 16 */ 17 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 18 19 /* 20 SuperLU include files 21 */ 22 EXTERN_C_BEGIN 23 #if defined(PETSC_USE_COMPLEX) 24 #include <slu_zdefs.h> 25 #else 26 #include <slu_ddefs.h> 27 #endif 28 #include <slu_util.h> 29 EXTERN_C_END 30 31 /* 32 This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 33 */ 34 typedef struct { 35 SuperMatrix A,L,U,B,X; 36 superlu_options_t options; 37 PetscInt *perm_c; /* column permutation vector */ 38 PetscInt *perm_r; /* row permutations from partial pivoting */ 39 PetscInt *etree; 40 PetscReal *R, *C; 41 char equed[1]; 42 PetscInt lwork; 43 void *work; 44 PetscReal rpg, rcond; 45 mem_usage_t mem_usage; 46 MatStructure flg; 47 SuperLUStat_t stat; 48 Mat A_dup; 49 PetscScalar *rhs_dup; 50 51 /* Flag to clean up (non-global) SuperLU objects during Destroy */ 52 PetscBool CleanUpSuperLU; 53 } Mat_SuperLU; 54 55 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 56 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *); 57 extern PetscErrorCode MatDestroy_SuperLU(Mat); 58 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 59 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 60 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 61 extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 62 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 63 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 64 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *); 65 66 /* 67 Utility function 68 */ 69 #undef __FUNCT__ 70 #define __FUNCT__ "MatFactorInfo_SuperLU" 71 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 72 { 73 Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 74 PetscErrorCode ierr; 75 superlu_options_t options; 76 77 PetscFunctionBegin; 78 /* check if matrix is superlu_dist type */ 79 if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 80 81 options = lu->options; 82 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 83 ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 84 ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 87 ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 88 ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 89 ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 90 ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 91 ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 92 ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 93 ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 94 if (A->factortype == MAT_FACTOR_ILU){ 95 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 97 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 98 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 99 ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 100 ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 101 } 102 PetscFunctionReturn(0); 103 } 104 105 /* 106 These are the methods provided to REPLACE the corresponding methods of the 107 SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 108 */ 109 #undef __FUNCT__ 110 #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 111 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 112 { 113 Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 114 Mat_SeqAIJ *aa; 115 PetscErrorCode ierr; 116 PetscInt sinfo; 117 PetscReal ferr, berr; 118 NCformat *Ustore; 119 SCformat *Lstore; 120 121 PetscFunctionBegin; 122 if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 123 lu->options.Fact = SamePattern; 124 /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 125 Destroy_SuperMatrix_Store(&lu->A); 126 if (lu->options.Equil){ 127 ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 128 } 129 if ( lu->lwork >= 0 ) { 130 Destroy_SuperNode_Matrix(&lu->L); 131 Destroy_CompCol_Matrix(&lu->U); 132 lu->options.Fact = SamePattern; 133 } 134 } 135 136 /* Create the SuperMatrix for lu->A=A^T: 137 Since SuperLU likes column-oriented matrices,we pass it the transpose, 138 and then solve A^T X = B in MatSolve(). */ 139 if (lu->options.Equil){ 140 aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 141 } else { 142 aa = (Mat_SeqAIJ*)(A)->data; 143 } 144 #if defined(PETSC_USE_COMPLEX) 145 zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 146 SLU_NC,SLU_Z,SLU_GE); 147 #else 148 dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i, 149 SLU_NC,SLU_D,SLU_GE); 150 #endif 151 152 /* Numerical factorization */ 153 lu->B.ncol = 0; /* Indicate not to solve the system */ 154 if (F->factortype == MAT_FACTOR_LU){ 155 #if defined(PETSC_USE_COMPLEX) 156 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 157 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 158 &lu->mem_usage, &lu->stat, &sinfo); 159 #else 160 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 161 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 162 &lu->mem_usage, &lu->stat, &sinfo); 163 #endif 164 } else if (F->factortype == MAT_FACTOR_ILU){ 165 /* Compute the incomplete factorization, condition number and pivot growth */ 166 #if defined(PETSC_USE_COMPLEX) 167 zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 168 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 169 &lu->mem_usage, &lu->stat, &sinfo); 170 #else 171 dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 172 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 173 &lu->mem_usage, &lu->stat, &sinfo); 174 #endif 175 } else { 176 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 177 } 178 if ( !sinfo || sinfo == lu->A.ncol+1 ) { 179 if ( lu->options.PivotGrowth ) 180 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 181 if ( lu->options.ConditionNumber ) 182 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 183 } else if ( sinfo > 0 ){ 184 if ( lu->lwork == -1 ) { 185 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 186 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 187 } else { /* sinfo < 0 */ 188 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 189 } 190 191 if ( lu->options.PrintStat ) { 192 ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 193 StatPrint(&lu->stat); 194 Lstore = (SCformat *) lu->L.Store; 195 Ustore = (NCformat *) lu->U.Store; 196 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 197 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 198 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 199 ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 200 lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 201 } 202 203 lu->flg = SAME_NONZERO_PATTERN; 204 F->ops->solve = MatSolve_SuperLU; 205 F->ops->solvetranspose = MatSolveTranspose_SuperLU; 206 F->ops->matsolve = MatMatSolve_SuperLU; 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "MatDestroy_SuperLU" 212 PetscErrorCode MatDestroy_SuperLU(Mat A) 213 { 214 PetscErrorCode ierr; 215 Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 216 217 PetscFunctionBegin; 218 if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 219 Destroy_SuperMatrix_Store(&lu->A); 220 Destroy_SuperMatrix_Store(&lu->B); 221 Destroy_SuperMatrix_Store(&lu->X); 222 StatFree(&lu->stat); 223 if (lu->lwork >= 0) { 224 Destroy_SuperNode_Matrix(&lu->L); 225 Destroy_CompCol_Matrix(&lu->U); 226 } 227 } 228 if (lu) { 229 ierr = PetscFree(lu->etree);CHKERRQ(ierr); 230 ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 231 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 232 ierr = PetscFree(lu->R);CHKERRQ(ierr); 233 ierr = PetscFree(lu->C);CHKERRQ(ierr); 234 ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 235 ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr); 236 } 237 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 238 239 /* clear composed functions */ 240 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 241 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr); 242 243 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 244 PetscFunctionReturn(0); 245 } 246 247 #undef __FUNCT__ 248 #define __FUNCT__ "MatView_SuperLU" 249 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 250 { 251 PetscErrorCode ierr; 252 PetscBool iascii; 253 PetscViewerFormat format; 254 255 PetscFunctionBegin; 256 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 257 if (iascii) { 258 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 259 if (format == PETSC_VIEWER_ASCII_INFO) { 260 ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 261 } 262 } 263 PetscFunctionReturn(0); 264 } 265 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "MatSolve_SuperLU_Private" 269 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 270 { 271 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 272 PetscScalar *barray,*xarray; 273 PetscErrorCode ierr; 274 PetscInt info,i,n=x->map->n; 275 PetscReal ferr,berr; 276 277 PetscFunctionBegin; 278 if ( lu->lwork == -1 ) { 279 PetscFunctionReturn(0); 280 } 281 282 lu->B.ncol = 1; /* Set the number of right-hand side */ 283 if (lu->options.Equil && !lu->rhs_dup){ 284 /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 285 ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr); 286 } 287 if (lu->options.Equil){ 288 /* Copy b into rsh_dup */ 289 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 290 ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 291 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 292 barray = lu->rhs_dup; 293 } else { 294 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 295 } 296 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 297 298 #if defined(PETSC_USE_COMPLEX) 299 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 300 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 301 #else 302 ((DNformat*)lu->B.Store)->nzval = barray; 303 ((DNformat*)lu->X.Store)->nzval = xarray; 304 #endif 305 306 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 307 if (A->factortype == MAT_FACTOR_LU){ 308 #if defined(PETSC_USE_COMPLEX) 309 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 310 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 311 &lu->mem_usage, &lu->stat, &info); 312 #else 313 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 314 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 315 &lu->mem_usage, &lu->stat, &info); 316 #endif 317 } else if (A->factortype == MAT_FACTOR_ILU){ 318 #if defined(PETSC_USE_COMPLEX) 319 zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 320 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 321 &lu->mem_usage, &lu->stat, &info); 322 #else 323 dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 324 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 325 &lu->mem_usage, &lu->stat, &info); 326 #endif 327 } else { 328 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 329 } 330 if (!lu->options.Equil){ 331 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 332 } 333 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 334 335 if ( !info || info == lu->A.ncol+1 ) { 336 if ( lu->options.IterRefine ) { 337 ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 338 ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 339 for (i = 0; i < 1; ++i) 340 ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 341 } 342 } else if ( info > 0 ){ 343 if ( lu->lwork == -1 ) { 344 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 345 } else { 346 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 347 } 348 } else if (info < 0){ 349 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 350 } 351 352 if ( lu->options.PrintStat ) { 353 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 354 StatPrint(&lu->stat); 355 } 356 PetscFunctionReturn(0); 357 } 358 359 #undef __FUNCT__ 360 #define __FUNCT__ "MatSolve_SuperLU" 361 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 362 { 363 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 lu->options.Trans = TRANS; 368 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372 #undef __FUNCT__ 373 #define __FUNCT__ "MatSolveTranspose_SuperLU" 374 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 375 { 376 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 377 PetscErrorCode ierr; 378 379 PetscFunctionBegin; 380 lu->options.Trans = NOTRANS; 381 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "MatMatSolve_SuperLU" 387 PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 388 { 389 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 390 PetscBool flg; 391 PetscErrorCode ierr; 392 393 PetscFunctionBegin; 394 ierr = PetscTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 395 if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 396 ierr = PetscTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 397 if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); lu->options.Trans = TRANS; 398 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 399 PetscFunctionReturn(0); 400 } 401 402 /* 403 Note the r permutation is ignored 404 */ 405 #undef __FUNCT__ 406 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 407 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 408 { 409 Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 410 411 PetscFunctionBegin; 412 lu->flg = DIFFERENT_NONZERO_PATTERN; 413 lu->CleanUpSuperLU = PETSC_TRUE; 414 F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 415 PetscFunctionReturn(0); 416 } 417 418 EXTERN_C_BEGIN 419 #undef __FUNCT__ 420 #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 421 PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 422 { 423 Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 424 425 PetscFunctionBegin; 426 lu->options.ILU_DropTol = dtol; 427 PetscFunctionReturn(0); 428 } 429 EXTERN_C_END 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "MatSuperluSetILUDropTol" 433 /*@ 434 MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 435 Logically Collective on Mat 436 437 Input Parameters: 438 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 439 - dtol - drop tolerance 440 441 Options Database: 442 . -mat_superlu_ilu_droptol <dtol> 443 444 Level: beginner 445 446 References: SuperLU Users' Guide 447 448 .seealso: MatGetFactor() 449 @*/ 450 PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 451 { 452 PetscErrorCode ierr; 453 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 456 PetscValidLogicalCollectiveInt(F,dtol,2); 457 ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460 461 EXTERN_C_BEGIN 462 #undef __FUNCT__ 463 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 464 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 465 { 466 PetscFunctionBegin; 467 *type = MATSOLVERSUPERLU; 468 PetscFunctionReturn(0); 469 } 470 EXTERN_C_END 471 472 473 /*MC 474 MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 475 via the external package SuperLU. 476 477 Use ./configure --download-superlu to have PETSc installed with SuperLU 478 479 Options Database Keys: 480 + -mat_superlu_equil: <FALSE> Equil (None) 481 . -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 482 . -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA 483 . -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None) 484 . -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None) 485 . -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None) 486 . -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None) 487 . -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag 488 . -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None) 489 . -mat_superlu_printstat: <FALSE> PrintStat (None) 490 . -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None) 491 . -mat_superlu_ilu_droptol <0>: ILU_DropTol (None) 492 . -mat_superlu_ilu_filltol <0>: ILU_FillTol (None) 493 . -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None) 494 . -mat_superlu_ilu_droprull <0>: ILU_DropRule (None) 495 . -mat_superlu_ilu_norm <0>: ILU_Norm (None) 496 - -mat_superlu_ilu_milu <0>: ILU_MILU (None) 497 498 Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 499 500 Level: beginner 501 502 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage 503 M*/ 504 505 EXTERN_C_BEGIN 506 #undef __FUNCT__ 507 #define __FUNCT__ "MatGetFactor_seqaij_superlu" 508 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 509 { 510 Mat B; 511 Mat_SuperLU *lu; 512 PetscErrorCode ierr; 513 PetscInt indx,m=A->rmap->n,n=A->cmap->n; 514 PetscBool flg; 515 const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 516 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 517 const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 518 519 PetscFunctionBegin; 520 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 521 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 522 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 523 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 524 525 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){ 526 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 527 B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 528 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 529 530 B->ops->destroy = MatDestroy_SuperLU; 531 B->ops->view = MatView_SuperLU; 532 B->factortype = ftype; 533 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 534 B->preallocated = PETSC_TRUE; 535 536 ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 537 538 if (ftype == MAT_FACTOR_LU){ 539 set_default_options(&lu->options); 540 /* Comments from SuperLU_4.0/SRC/dgssvx.c: 541 "Whether or not the system will be equilibrated depends on the 542 scaling of the matrix A, but if equilibration is used, A is 543 overwritten by diag(R)*A*diag(C) and B by diag(R)*B 544 (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 545 We set 'options.Equil = NO' as default because additional space is needed for it. 546 */ 547 lu->options.Equil = NO; 548 } else if (ftype == MAT_FACTOR_ILU){ 549 /* Set the default input options of ilu: */ 550 #if 0 551 options.Fact = DOFACT; 552 options.Equil = YES; /* must be YES for ilu - don't know why */ 553 options.ColPerm = COLAMD; 554 options.DiagPivotThresh = 0.1; /* different from complete LU */ 555 options.Trans = NOTRANS; 556 options.IterRefine = NOREFINE; 557 options.SymmetricMode = NO; 558 options.PivotGrowth = NO; 559 options.ConditionNumber = NO; 560 options.PrintStat = YES; 561 options.RowPerm = LargeDiag; 562 options.ILU_DropTol = 1e-4; 563 options.ILU_FillTol = 1e-2; 564 options.ILU_FillFactor = 10.0; 565 options.ILU_DropRule = DROP_BASIC | DROP_AREA; 566 options.ILU_Norm = INF_NORM; 567 options.ILU_MILU = SMILU_2; 568 #endif 569 ilu_set_default_options(&lu->options); 570 /* there is a bug with options.RowPerm=LargeDiag causing src/ksp/ksp/examples/tutorials/runex52_superlu crashes 571 See email communication Betwen Hong and Sharry on Feb. Tue, Feb 22, 2011 */ 572 lu->options.RowPerm = NOROWPERM; 573 } 574 lu->options.PrintStat = NO; 575 576 /* Initialize the statistics variables. */ 577 StatInit(&lu->stat); 578 lu->lwork = 0; /* allocate space internally by system malloc */ 579 580 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 581 ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 582 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 583 if (flg) {lu->options.ColPerm = (colperm_t)indx;} 584 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 585 if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 586 ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 587 if (flg) lu->options.SymmetricMode = YES; 588 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 589 ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 590 if (flg) lu->options.PivotGrowth = YES; 591 ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 592 if (flg) lu->options.ConditionNumber = YES; 593 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 594 if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 595 ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 596 if (flg) lu->options.ReplaceTinyPivot = YES; 597 ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 598 if (flg) lu->options.PrintStat = YES; 599 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 600 if (lu->lwork > 0 ){ 601 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 602 } else if (lu->lwork != 0 && lu->lwork != -1){ 603 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 604 lu->lwork = 0; 605 } 606 /* ilu options */ 607 ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr); 608 ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr); 609 ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr); 610 ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 611 ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 612 if (flg){ 613 lu->options.ILU_Norm = (norm_t)indx; 614 } 615 ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 616 if (flg){ 617 lu->options.ILU_MILU = (milu_t)indx; 618 } 619 PetscOptionsEnd(); 620 if (lu->options.Equil == YES) { 621 /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 622 ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 623 } 624 625 /* Allocate spaces (notice sizes are for the transpose) */ 626 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 627 ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 628 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 629 ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 630 ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 631 632 /* create rhs and solution x without allocate space for .Store */ 633 #if defined(PETSC_USE_COMPLEX) 634 zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 635 zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 636 #else 637 dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 638 dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 639 #endif 640 641 #ifdef SUPERLU2 642 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 643 #endif 644 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 645 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 646 B->spptr = lu; 647 *F = B; 648 PetscFunctionReturn(0); 649 } 650 EXTERN_C_END 651 652