1 #define PETSCMAT_DLL 2 3 /* -------------------------------------------------------------------- 4 5 This file implements a subclass of the SeqAIJ matrix class that uses 6 the SuperLU sparse solver. You can use this as a starting point for 7 implementing your own subclass of a PETSc matrix class. 8 9 This demonstrates a way to make an implementation inheritence of a PETSc 10 matrix type. This means constructing a new matrix type (SuperLU) by changing some 11 of the methods of the previous type (SeqAIJ), adding additional data, and possibly 12 additional method. (See any book on object oriented programming). 13 */ 14 15 /* 16 Defines the data structure for the base matrix type (SeqAIJ) 17 */ 18 #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 19 20 /* 21 SuperLU include files 22 */ 23 EXTERN_C_BEGIN 24 #if defined(PETSC_USE_COMPLEX) 25 #include "slu_zdefs.h" 26 #else 27 #include "slu_ddefs.h" 28 #endif 29 #include "slu_util.h" 30 EXTERN_C_END 31 32 /* 33 This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 34 */ 35 typedef struct { 36 SuperMatrix A,L,U,B,X; 37 superlu_options_t options; 38 PetscInt *perm_c; /* column permutation vector */ 39 PetscInt *perm_r; /* row permutations from partial pivoting */ 40 PetscInt *etree; 41 PetscReal *R, *C; 42 char equed[1]; 43 PetscInt lwork; 44 void *work; 45 PetscReal rpg, rcond; 46 mem_usage_t mem_usage; 47 MatStructure flg; 48 SuperLUStat_t stat; 49 Mat A_dup; 50 PetscScalar *rhs_dup; 51 52 /* Flag to clean up (non-global) SuperLU objects during Destroy */ 53 PetscBool CleanUpSuperLU; 54 } Mat_SuperLU; 55 56 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 57 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *); 58 extern PetscErrorCode MatDestroy_SuperLU(Mat); 59 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 60 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 61 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 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 PetscFunctionReturn(0); 207 } 208 209 #undef __FUNCT__ 210 #define __FUNCT__ "MatDestroy_SuperLU" 211 PetscErrorCode MatDestroy_SuperLU(Mat A) 212 { 213 PetscErrorCode ierr; 214 Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 215 216 PetscFunctionBegin; 217 if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 218 Destroy_SuperMatrix_Store(&lu->A); 219 Destroy_SuperMatrix_Store(&lu->B); 220 Destroy_SuperMatrix_Store(&lu->X); 221 StatFree(&lu->stat); 222 if ( lu->lwork >= 0 ) { 223 Destroy_SuperNode_Matrix(&lu->L); 224 Destroy_CompCol_Matrix(&lu->U); 225 } 226 227 } 228 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 235 /* clear composed functions */ 236 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 237 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr); 238 239 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 240 if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);} 241 if (lu->rhs_dup){ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);} 242 PetscFunctionReturn(0); 243 } 244 245 #undef __FUNCT__ 246 #define __FUNCT__ "MatView_SuperLU" 247 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 248 { 249 PetscErrorCode ierr; 250 PetscBool iascii; 251 PetscViewerFormat format; 252 253 PetscFunctionBegin; 254 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 255 if (iascii) { 256 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 257 if (format == PETSC_VIEWER_ASCII_INFO) { 258 ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 259 } 260 } 261 PetscFunctionReturn(0); 262 } 263 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "MatSolve_SuperLU_Private" 267 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 268 { 269 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 270 PetscScalar *barray,*xarray; 271 PetscErrorCode ierr; 272 PetscInt info,i,n=x->map->n; 273 PetscReal ferr,berr; 274 275 PetscFunctionBegin; 276 if ( lu->lwork == -1 ) { 277 PetscFunctionReturn(0); 278 } 279 280 lu->B.ncol = 1; /* Set the number of right-hand side */ 281 if (lu->options.Equil && !lu->rhs_dup){ 282 /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 283 ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr); 284 } 285 if (lu->options.Equil){ 286 /* Copy b into rsh_dup */ 287 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 288 ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 289 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 290 barray = lu->rhs_dup; 291 } else { 292 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 293 } 294 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 295 296 #if defined(PETSC_USE_COMPLEX) 297 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 298 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 299 #else 300 ((DNformat*)lu->B.Store)->nzval = barray; 301 ((DNformat*)lu->X.Store)->nzval = xarray; 302 #endif 303 304 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 305 if (A->factortype == MAT_FACTOR_LU){ 306 #if defined(PETSC_USE_COMPLEX) 307 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 308 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 309 &lu->mem_usage, &lu->stat, &info); 310 #else 311 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 312 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 313 &lu->mem_usage, &lu->stat, &info); 314 #endif 315 } else if (A->factortype == MAT_FACTOR_ILU){ 316 #if defined(PETSC_USE_COMPLEX) 317 zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 318 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 319 &lu->mem_usage, &lu->stat, &info); 320 #else 321 dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 322 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 323 &lu->mem_usage, &lu->stat, &info); 324 #endif 325 } else { 326 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 327 } 328 if (!lu->options.Equil){ 329 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 330 } 331 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 332 333 if ( !info || info == lu->A.ncol+1 ) { 334 if ( lu->options.IterRefine ) { 335 ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 336 ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 337 for (i = 0; i < 1; ++i) 338 ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 339 } 340 } else if ( info > 0 ){ 341 if ( lu->lwork == -1 ) { 342 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 343 } else { 344 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 345 } 346 } else if (info < 0){ 347 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 348 } 349 350 if ( lu->options.PrintStat ) { 351 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 352 StatPrint(&lu->stat); 353 } 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "MatSolve_SuperLU" 359 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 360 { 361 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 lu->options.Trans = TRANS; 366 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 367 PetscFunctionReturn(0); 368 } 369 370 #undef __FUNCT__ 371 #define __FUNCT__ "MatSolveTranspose_SuperLU" 372 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 373 { 374 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 375 PetscErrorCode ierr; 376 377 PetscFunctionBegin; 378 lu->options.Trans = NOTRANS; 379 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 380 PetscFunctionReturn(0); 381 } 382 383 /* 384 Note the r permutation is ignored 385 */ 386 #undef __FUNCT__ 387 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 388 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 389 { 390 Mat_SuperLU *lu = (Mat_SuperLU*)((F)->spptr); 391 392 PetscFunctionBegin; 393 lu->flg = DIFFERENT_NONZERO_PATTERN; 394 lu->CleanUpSuperLU = PETSC_TRUE; 395 (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 396 PetscFunctionReturn(0); 397 } 398 399 EXTERN_C_BEGIN 400 #undef __FUNCT__ 401 #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 402 PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 403 { 404 Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 405 406 PetscFunctionBegin; 407 lu->options.ILU_DropTol = dtol; 408 PetscFunctionReturn(0); 409 } 410 EXTERN_C_END 411 412 #undef __FUNCT__ 413 #define __FUNCT__ "MatSuperluSetILUDropTol" 414 /*@ 415 MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 416 Logically Collective on Mat 417 418 Input Parameters: 419 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 420 - dtol - drop tolerance 421 422 Options Database: 423 . -mat_superlu_ilu_droptol <dtol> 424 425 Level: beginner 426 427 References: SuperLU Users' Guide 428 429 .seealso: MatGetFactor() 430 @*/ 431 PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 432 { 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 437 PetscValidLogicalCollectiveInt(F,dtol,2); 438 ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 EXTERN_C_BEGIN 443 #undef __FUNCT__ 444 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 445 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 446 { 447 PetscFunctionBegin; 448 *type = MATSOLVERSUPERLU; 449 PetscFunctionReturn(0); 450 } 451 EXTERN_C_END 452 453 454 /*MC 455 MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 456 via the external package SuperLU. 457 458 Use ./configure --download-superlu to have PETSc installed with SuperLU 459 460 Options Database Keys: 461 + -mat_superlu_equil: <FALSE> Equil (None) 462 . -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 463 . -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA 464 . -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None) 465 . -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None) 466 . -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None) 467 . -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None) 468 . -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag 469 . -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None) 470 . -mat_superlu_printstat: <FALSE> PrintStat (None) 471 . -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None) 472 . -mat_superlu_ilu_droptol <0>: ILU_DropTol (None) 473 . -mat_superlu_ilu_filltol <0>: ILU_FillTol (None) 474 . -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None) 475 . -mat_superlu_ilu_droprull <0>: ILU_DropRule (None) 476 . -mat_superlu_ilu_norm <0>: ILU_Norm (None) 477 - -mat_superlu_ilu_milu <0>: ILU_MILU (None) 478 479 Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 480 481 Level: beginner 482 483 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage 484 M*/ 485 486 EXTERN_C_BEGIN 487 #undef __FUNCT__ 488 #define __FUNCT__ "MatGetFactor_seqaij_superlu" 489 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 490 { 491 Mat B; 492 Mat_SuperLU *lu; 493 PetscErrorCode ierr; 494 PetscInt indx,m=A->rmap->n,n=A->cmap->n; 495 PetscBool flg; 496 const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 497 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 498 const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 499 500 PetscFunctionBegin; 501 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 502 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 503 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 504 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 505 506 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){ 507 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 508 B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 509 } else { 510 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 511 } 512 513 B->ops->destroy = MatDestroy_SuperLU; 514 B->ops->view = MatView_SuperLU; 515 B->factortype = ftype; 516 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 517 B->preallocated = PETSC_TRUE; 518 519 ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 520 521 if (ftype == MAT_FACTOR_LU){ 522 set_default_options(&lu->options); 523 /* Comments from SuperLU_4.0/SRC/dgssvx.c: 524 "Whether or not the system will be equilibrated depends on the 525 scaling of the matrix A, but if equilibration is used, A is 526 overwritten by diag(R)*A*diag(C) and B by diag(R)*B 527 (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 528 We set 'options.Equil = NO' as default because additional space is needed for it. 529 */ 530 lu->options.Equil = NO; 531 } else if (ftype == MAT_FACTOR_ILU){ 532 /* Set the default input options of ilu: 533 options.Fact = DOFACT; 534 options.Equil = YES; // must be YES for ilu - don't know why 535 options.ColPerm = COLAMD; 536 options.DiagPivotThresh = 0.1; //different from complete LU 537 options.Trans = NOTRANS; 538 options.IterRefine = NOREFINE; 539 options.SymmetricMode = NO; 540 options.PivotGrowth = NO; 541 options.ConditionNumber = NO; 542 options.PrintStat = YES; 543 options.RowPerm = LargeDiag; 544 options.ILU_DropTol = 1e-4; 545 options.ILU_FillTol = 1e-2; 546 options.ILU_FillFactor = 10.0; 547 options.ILU_DropRule = DROP_BASIC | DROP_AREA; 548 options.ILU_Norm = INF_NORM; 549 options.ILU_MILU = SMILU_2; 550 */ 551 ilu_set_default_options(&lu->options); 552 } 553 lu->options.PrintStat = NO; 554 555 /* Initialize the statistics variables. */ 556 StatInit(&lu->stat); 557 lu->lwork = 0; /* allocate space internally by system malloc */ 558 559 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 560 ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 561 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 562 if (flg) {lu->options.ColPerm = (colperm_t)indx;} 563 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 564 if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 565 ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 566 if (flg) lu->options.SymmetricMode = YES; 567 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 568 ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 569 if (flg) lu->options.PivotGrowth = YES; 570 ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 571 if (flg) lu->options.ConditionNumber = YES; 572 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 573 if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 574 ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 575 if (flg) lu->options.ReplaceTinyPivot = YES; 576 ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 577 if (flg) lu->options.PrintStat = YES; 578 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 579 if (lu->lwork > 0 ){ 580 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 581 } else if (lu->lwork != 0 && lu->lwork != -1){ 582 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 583 lu->lwork = 0; 584 } 585 /* ilu options */ 586 ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr); 587 ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr); 588 ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr); 589 ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 590 ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 591 if (flg){ 592 lu->options.ILU_Norm = (norm_t)indx; 593 } 594 ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 595 if (flg){ 596 lu->options.ILU_MILU = (milu_t)indx; 597 } 598 PetscOptionsEnd(); 599 if (lu->options.Equil == YES) { 600 /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 601 ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 602 } 603 604 /* Allocate spaces (notice sizes are for the transpose) */ 605 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 606 ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 607 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 608 ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 609 ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 610 611 /* create rhs and solution x without allocate space for .Store */ 612 #if defined(PETSC_USE_COMPLEX) 613 zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 614 zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 615 #else 616 dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 617 dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 618 #endif 619 620 #ifdef SUPERLU2 621 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 622 #endif 623 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 624 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 625 B->spptr = lu; 626 *F = B; 627 PetscFunctionReturn(0); 628 } 629 EXTERN_C_END 630 631