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