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