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