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" 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 50 /* Flag to clean up (non-global) SuperLU objects during Destroy */ 51 PetscTruth CleanUpSuperLU; 52 } Mat_SuperLU; 53 54 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 55 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *); 56 extern PetscErrorCode MatDestroy_SuperLU(Mat); 57 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 58 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 59 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 60 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 61 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 62 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *); 63 64 /* 65 Utility function 66 */ 67 #undef __FUNCT__ 68 #define __FUNCT__ "MatFactorInfo_SuperLU" 69 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 70 { 71 Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 72 PetscErrorCode ierr; 73 superlu_options_t options; 74 75 PetscFunctionBegin; 76 /* check if matrix is superlu_dist type */ 77 if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 78 79 options = lu->options; 80 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 81 ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 82 ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 83 ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 84 ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 87 ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 88 ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 89 ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 90 ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 91 ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 92 if (A->factor == MAT_FACTOR_ILU){ 93 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 94 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 95 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 97 ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 98 ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 99 } 100 PetscFunctionReturn(0); 101 } 102 103 /* 104 These are the methods provided to REPLACE the corresponding methods of the 105 SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 106 */ 107 #undef __FUNCT__ 108 #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 109 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 110 { 111 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 112 Mat_SuperLU *lu = (Mat_SuperLU*)(F)->spptr; 113 PetscErrorCode ierr; 114 PetscInt sinfo; 115 PetscReal ferr, berr; 116 NCformat *Ustore; 117 SCformat *Lstore; 118 119 PetscFunctionBegin; 120 if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 121 lu->options.Fact = SamePattern; 122 /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 123 Destroy_SuperMatrix_Store(&lu->A); 124 if ( lu->lwork >= 0 ) { 125 Destroy_SuperNode_Matrix(&lu->L); 126 Destroy_CompCol_Matrix(&lu->U); 127 lu->options.Fact = SamePattern; 128 } 129 } 130 131 /* Create the SuperMatrix for lu->A=A^T: 132 Since SuperLU likes column-oriented matrices,we pass it the transpose, 133 and then solve A^T X = B in MatSolve(). */ 134 #if defined(PETSC_USE_COMPLEX) 135 zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 136 SLU_NC,SLU_Z,SLU_GE); 137 #else 138 dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i, 139 SLU_NC,SLU_D,SLU_GE); 140 #endif 141 142 /* Numerical factorization */ 143 lu->B.ncol = 0; /* Indicate not to solve the system */ 144 if (F->factor == MAT_FACTOR_LU){ 145 #if defined(PETSC_USE_COMPLEX) 146 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 147 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 148 &lu->mem_usage, &lu->stat, &sinfo); 149 #else 150 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 151 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 152 &lu->mem_usage, &lu->stat, &sinfo); 153 #endif 154 } else if (F->factor == MAT_FACTOR_ILU){ 155 /* Compute the incomplete factorization, condition number and pivot growth */ 156 #if defined(PETSC_USE_COMPLEX) 157 zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 158 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 159 &lu->mem_usage, &lu->stat, &sinfo); 160 #else 161 dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 162 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 163 &lu->mem_usage, &lu->stat, &sinfo); 164 #endif 165 } else { 166 SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 167 } 168 if ( !sinfo || sinfo == lu->A.ncol+1 ) { 169 if ( lu->options.PivotGrowth ) 170 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 171 if ( lu->options.ConditionNumber ) 172 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 173 } else if ( sinfo > 0 ){ 174 if ( lu->lwork == -1 ) { 175 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 176 } else { 177 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 178 } 179 } else { /* sinfo < 0 */ 180 SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 181 } 182 183 if ( lu->options.PrintStat ) { 184 ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 185 StatPrint(&lu->stat); 186 Lstore = (SCformat *) lu->L.Store; 187 Ustore = (NCformat *) lu->U.Store; 188 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 189 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 190 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 191 ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 192 lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 193 } 194 195 lu->flg = SAME_NONZERO_PATTERN; 196 (F)->ops->solve = MatSolve_SuperLU; 197 (F)->ops->solvetranspose = MatSolveTranspose_SuperLU; 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "MatDestroy_SuperLU" 203 PetscErrorCode MatDestroy_SuperLU(Mat A) 204 { 205 PetscErrorCode ierr; 206 Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 207 208 PetscFunctionBegin; 209 if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 210 Destroy_SuperMatrix_Store(&lu->A); 211 Destroy_SuperMatrix_Store(&lu->B); 212 Destroy_SuperMatrix_Store(&lu->X); 213 StatFree(&lu->stat); 214 215 ierr = PetscFree(lu->etree);CHKERRQ(ierr); 216 ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 217 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 218 ierr = PetscFree(lu->R);CHKERRQ(ierr); 219 ierr = PetscFree(lu->C);CHKERRQ(ierr); 220 if ( lu->lwork >= 0 ) { 221 Destroy_SuperNode_Matrix(&lu->L); 222 Destroy_CompCol_Matrix(&lu->U); 223 } 224 } 225 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 226 PetscFunctionReturn(0); 227 } 228 229 #undef __FUNCT__ 230 #define __FUNCT__ "MatView_SuperLU" 231 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 232 { 233 PetscErrorCode ierr; 234 PetscTruth iascii; 235 PetscViewerFormat format; 236 237 PetscFunctionBegin; 238 ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr); 239 240 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 241 if (iascii) { 242 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 243 if (format == PETSC_VIEWER_ASCII_INFO) { 244 ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 245 } 246 } 247 PetscFunctionReturn(0); 248 } 249 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "MatSolve_SuperLU_Private" 253 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 254 { 255 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 256 PetscScalar *barray,*xarray; 257 PetscErrorCode ierr; 258 PetscInt info,i; 259 PetscReal ferr,berr; 260 261 PetscFunctionBegin; 262 if ( lu->lwork == -1 ) { 263 PetscFunctionReturn(0); 264 } 265 lu->B.ncol = 1; /* Set the number of right-hand side */ 266 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 267 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 268 269 #if defined(PETSC_USE_COMPLEX) 270 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 271 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 272 #else 273 ((DNformat*)lu->B.Store)->nzval = barray; 274 ((DNformat*)lu->X.Store)->nzval = xarray; 275 #endif 276 277 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 278 if (A->factor == MAT_FACTOR_LU){ 279 #if defined(PETSC_USE_COMPLEX) 280 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 281 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 282 &lu->mem_usage, &lu->stat, &info); 283 #else 284 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 285 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 286 &lu->mem_usage, &lu->stat, &info); 287 #endif 288 } else if (A->factor == MAT_FACTOR_ILU){ 289 #if defined(PETSC_USE_COMPLEX) 290 zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 291 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 292 &lu->mem_usage, &lu->stat, &info); 293 #else 294 dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 295 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 296 &lu->mem_usage, &lu->stat, &info); 297 #endif 298 } else { 299 SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 300 } 301 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 302 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 303 304 if ( !info || info == lu->A.ncol+1 ) { 305 if ( lu->options.IterRefine ) { 306 ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 307 ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 308 for (i = 0; i < 1; ++i) 309 ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 310 } 311 } else if ( info > 0 ){ 312 if ( lu->lwork == -1 ) { 313 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 314 } else { 315 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 316 } 317 } else if (info < 0){ 318 SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 319 } 320 321 if ( lu->options.PrintStat ) { 322 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 323 StatPrint(&lu->stat); 324 } 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "MatSolve_SuperLU" 330 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 331 { 332 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 333 PetscErrorCode ierr; 334 335 PetscFunctionBegin; 336 lu->options.Trans = TRANS; 337 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "MatSolveTranspose_SuperLU" 343 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 344 { 345 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 346 PetscErrorCode ierr; 347 348 PetscFunctionBegin; 349 lu->options.Trans = NOTRANS; 350 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 351 PetscFunctionReturn(0); 352 } 353 354 /* 355 Note the r permutation is ignored 356 */ 357 #undef __FUNCT__ 358 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 359 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 360 { 361 Mat_SuperLU *lu = (Mat_SuperLU*)((F)->spptr); 362 363 PetscFunctionBegin; 364 lu->flg = DIFFERENT_NONZERO_PATTERN; 365 lu->CleanUpSuperLU = PETSC_TRUE; 366 (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 367 PetscFunctionReturn(0); 368 } 369 370 EXTERN_C_BEGIN 371 #undef __FUNCT__ 372 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 373 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 374 { 375 PetscFunctionBegin; 376 *type = MAT_SOLVER_SUPERLU; 377 PetscFunctionReturn(0); 378 } 379 EXTERN_C_END 380 381 382 /*MC 383 MAT_SOLVER_SUPERLU = "superlu" - A solver package roviding direct solvers (LU) for sequential matrices 384 via the external package SuperLU. 385 386 Use config/configure.py --download-superlu to have PETSc installed with SuperLU 387 388 Options Database Keys: 389 + -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 390 1: MMD applied to A'*A, 391 2: MMD applied to A'+A, 392 3: COLAMD, approximate minimum degree column ordering 393 . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve 394 choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE 395 - -mat_superlu_printstat - print SuperLU statistics about the factorization 396 397 Notes: Do not confuse this with MAT_SOLVER_SUPERLU_DIST which is for parallel sparse solves 398 399 Level: beginner 400 401 .seealso: PCLU, MAT_SOLVER_SUPERLU_DIST, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage 402 M*/ 403 404 EXTERN_C_BEGIN 405 #undef __FUNCT__ 406 #define __FUNCT__ "MatGetFactor_seqaij_superlu" 407 PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 408 { 409 Mat B; 410 Mat_SuperLU *lu; 411 PetscErrorCode ierr; 412 PetscInt indx,m=A->rmap->n,n=A->cmap->n; 413 PetscTruth flg; 414 const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 415 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 416 const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 417 418 PetscFunctionBegin; 419 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 420 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 421 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 422 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 423 424 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){ 425 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 426 B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 427 } else { 428 SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 429 } 430 431 B->ops->destroy = MatDestroy_SuperLU; 432 B->ops->view = MatView_SuperLU; 433 B->factor = ftype; 434 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 435 B->preallocated = PETSC_TRUE; 436 437 ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 438 if (ftype == MAT_FACTOR_LU){ 439 set_default_options(&lu->options); 440 } else if (ftype == MAT_FACTOR_ILU){ 441 /* Set the default input options of ilu: 442 options.Fact = DOFACT; 443 options.Equil = YES; 444 options.ColPerm = COLAMD; 445 options.DiagPivotThresh = 0.1; //different from complete LU 446 options.Trans = NOTRANS; 447 options.IterRefine = NOREFINE; 448 options.SymmetricMode = NO; 449 options.PivotGrowth = NO; 450 options.ConditionNumber = NO; 451 options.PrintStat = YES; 452 options.RowPerm = LargeDiag; 453 options.ILU_DropTol = 1e-4; 454 options.ILU_FillTol = 1e-2; 455 options.ILU_FillFactor = 10.0; 456 options.ILU_DropRule = DROP_BASIC | DROP_AREA; 457 options.ILU_Norm = INF_NORM; 458 options.ILU_MILU = SMILU_2; 459 */ 460 461 ilu_set_default_options(&lu->options); 462 } 463 /* equilibration causes error in solve()(ref. [petsc-maint #42782] SuperLU troubles) 464 thus not supported here. See dgssvx.c for possible reason. */ 465 lu->options.Equil = NO; 466 lu->options.PrintStat = NO; 467 468 /* Initialize the statistics variables. */ 469 StatInit(&lu->stat); 470 lu->lwork = 0; /* allocate space internally by system malloc */ 471 472 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 473 ierr = PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 474 if (flg) lu->options.Equil = YES; 475 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 476 if (flg) {lu->options.ColPerm = (colperm_t)indx;} 477 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 478 if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 479 ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 480 if (flg) lu->options.SymmetricMode = YES; 481 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 482 ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 483 if (flg) lu->options.PivotGrowth = YES; 484 ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 485 if (flg) lu->options.ConditionNumber = YES; 486 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 487 if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 488 ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 489 if (flg) lu->options.ReplaceTinyPivot = YES; 490 ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 491 if (flg) lu->options.PrintStat = YES; 492 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 493 if (lu->lwork > 0 ){ 494 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 495 } else if (lu->lwork != 0 && lu->lwork != -1){ 496 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 497 lu->lwork = 0; 498 } 499 /* ilu options */ 500 ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr); 501 ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr); 502 ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr); 503 ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 504 ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 505 if (flg){ 506 lu->options.ILU_Norm = (norm_t)indx; 507 } 508 ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 509 if (flg){ 510 lu->options.ILU_MILU = (milu_t)indx; 511 } 512 PetscOptionsEnd(); 513 514 /* Allocate spaces (notice sizes are for the transpose) */ 515 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 516 ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 517 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 518 ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 519 ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 520 521 /* create rhs and solution x without allocate space for .Store */ 522 #if defined(PETSC_USE_COMPLEX) 523 zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 524 zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 525 #else 526 dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 527 dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 528 #endif 529 530 #ifdef SUPERLU2 531 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 532 #endif 533 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 534 B->spptr = lu; 535 *F = B; 536 PetscFunctionReturn(0); 537 } 538 EXTERN_C_END 539