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