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