1 2 /* -------------------------------------------------------------------- 3 4 This file implements a subclass of the SeqAIJ matrix class that uses 5 the SuperLU sparse solver. 6 */ 7 8 /* 9 Defines the data structure for the base matrix type (SeqAIJ) 10 */ 11 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 12 13 /* 14 SuperLU include files 15 */ 16 EXTERN_C_BEGIN 17 #if defined(PETSC_USE_COMPLEX) 18 #if defined(PETSC_USE_REAL_SINGLE) 19 #include <slu_cdefs.h> 20 #else 21 #include <slu_zdefs.h> 22 #endif 23 #else 24 #if defined(PETSC_USE_REAL_SINGLE) 25 #include <slu_sdefs.h> 26 #else 27 #include <slu_ddefs.h> 28 #endif 29 #endif 30 #include <slu_util.h> 31 EXTERN_C_END 32 33 /* 34 This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 35 */ 36 typedef struct { 37 SuperMatrix A,L,U,B,X; 38 superlu_options_t options; 39 PetscInt *perm_c; /* column permutation vector */ 40 PetscInt *perm_r; /* row permutations from partial pivoting */ 41 PetscInt *etree; 42 PetscReal *R, *C; 43 char equed[1]; 44 PetscInt lwork; 45 void *work; 46 PetscReal rpg, rcond; 47 mem_usage_t mem_usage; 48 MatStructure flg; 49 SuperLUStat_t stat; 50 Mat A_dup; 51 PetscScalar *rhs_dup; 52 GlobalLU_t Glu; 53 54 /* Flag to clean up (non-global) SuperLU objects during Destroy */ 55 PetscBool CleanUpSuperLU; 56 } Mat_SuperLU; 57 58 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 59 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo*); 60 extern PetscErrorCode MatDestroy_SuperLU(Mat); 61 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 62 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 63 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 64 extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 65 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 66 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 67 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat*); 68 69 /* 70 Utility function 71 */ 72 #undef __FUNCT__ 73 #define __FUNCT__ "MatFactorInfo_SuperLU" 74 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 75 { 76 Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 77 PetscErrorCode ierr; 78 superlu_options_t options; 79 80 PetscFunctionBegin; 81 /* check if matrix is superlu_dist type */ 82 if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 83 84 options = lu->options; 85 86 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 87 ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");CHKERRQ(ierr); 88 ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 89 ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 90 ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");CHKERRQ(ierr); 91 ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 92 ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");CHKERRQ(ierr); 93 ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");CHKERRQ(ierr); 94 ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 95 ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");CHKERRQ(ierr); 97 ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 98 if (A->factortype == MAT_FACTOR_ILU) { 99 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 100 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 101 ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 102 ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 103 ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 104 ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 105 } 106 PetscFunctionReturn(0); 107 } 108 109 /* 110 These are the methods provided to REPLACE the corresponding methods of the 111 SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 112 */ 113 #undef __FUNCT__ 114 #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 115 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 116 { 117 Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 118 Mat_SeqAIJ *aa; 119 PetscErrorCode ierr; 120 PetscInt sinfo; 121 PetscReal ferr, berr; 122 NCformat *Ustore; 123 SCformat *Lstore; 124 125 PetscFunctionBegin; 126 if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */ 127 lu->options.Fact = SamePattern; 128 /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 129 Destroy_SuperMatrix_Store(&lu->A); 130 if (lu->options.Equil) { 131 ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 132 } 133 if (lu->lwork >= 0) { 134 PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 135 PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 136 lu->options.Fact = SamePattern; 137 } 138 } 139 140 /* Create the SuperMatrix for lu->A=A^T: 141 Since SuperLU likes column-oriented matrices,we pass it the transpose, 142 and then solve A^T X = B in MatSolve(). */ 143 if (lu->options.Equil) { 144 aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 145 } else { 146 aa = (Mat_SeqAIJ*)(A)->data; 147 } 148 #if defined(PETSC_USE_COMPLEX) 149 #if defined(PETSC_USE_REAL_SINGLE) 150 PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE)); 151 #else 152 PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE)); 153 #endif 154 #else 155 #if defined(PETSC_USE_REAL_SINGLE) 156 PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE)); 157 #else 158 PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE)); 159 #endif 160 #endif 161 162 /* Numerical factorization */ 163 lu->B.ncol = 0; /* Indicate not to solve the system */ 164 if (F->factortype == MAT_FACTOR_LU) { 165 #if defined(PETSC_USE_COMPLEX) 166 #if defined(PETSC_USE_REAL_SINGLE) 167 PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 168 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 169 &lu->mem_usage, &lu->stat, &sinfo)); 170 #else 171 PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 172 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 173 &lu->mem_usage, &lu->stat, &sinfo)); 174 #endif 175 #else 176 #if defined(PETSC_USE_REAL_SINGLE) 177 PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 178 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 179 &lu->mem_usage, &lu->stat, &sinfo)); 180 #else 181 PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 182 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 183 &lu->Glu,&lu->mem_usage, &lu->stat, &sinfo)); 184 #endif 185 #endif 186 } else if (F->factortype == MAT_FACTOR_ILU) { 187 /* Compute the incomplete factorization, condition number and pivot growth */ 188 #if defined(PETSC_USE_COMPLEX) 189 #if defined(PETSC_USE_REAL_SINGLE) 190 PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 191 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 192 &lu->mem_usage, &lu->stat, &sinfo)); 193 #else 194 PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 195 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 196 &lu->mem_usage, &lu->stat, &sinfo)); 197 #endif 198 #else 199 #if defined(PETSC_USE_REAL_SINGLE) 200 PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 201 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 202 &lu->mem_usage, &lu->stat, &sinfo)); 203 #else 204 PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 205 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 206 &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 207 #endif 208 #endif 209 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 210 if (!sinfo || sinfo == lu->A.ncol+1) { 211 if (lu->options.PivotGrowth) { 212 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 213 } 214 if (lu->options.ConditionNumber) { 215 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 216 } 217 } else if (sinfo > 0) { 218 if (lu->lwork == -1) { 219 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 220 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 221 } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 222 223 if (lu->options.PrintStat) { 224 ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 225 PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 226 Lstore = (SCformat*) lu->L.Store; 227 Ustore = (NCformat*) lu->U.Store; 228 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 229 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 230 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 231 ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 232 lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 233 } 234 235 lu->flg = SAME_NONZERO_PATTERN; 236 F->ops->solve = MatSolve_SuperLU; 237 F->ops->solvetranspose = MatSolveTranspose_SuperLU; 238 F->ops->matsolve = NULL; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "MatGetDiagonal_SuperLU" 244 PetscErrorCode MatGetDiagonal_SuperLU(Mat A,Vec v) 245 { 246 PetscFunctionBegin; 247 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: SuperLU factor"); 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "MatDestroy_SuperLU" 253 PetscErrorCode MatDestroy_SuperLU(Mat A) 254 { 255 PetscErrorCode ierr; 256 Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 257 258 PetscFunctionBegin; 259 if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 260 PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A)); 261 PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B)); 262 PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X)); 263 PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat)); 264 if (lu->lwork >= 0) { 265 PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 266 PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 267 } 268 } 269 if (lu) { 270 ierr = PetscFree(lu->etree);CHKERRQ(ierr); 271 ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 272 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 273 ierr = PetscFree(lu->R);CHKERRQ(ierr); 274 ierr = PetscFree(lu->C);CHKERRQ(ierr); 275 ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 276 ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr); 277 } 278 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 279 280 /* clear composed functions */ 281 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 282 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);CHKERRQ(ierr); 283 284 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "MatView_SuperLU" 290 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 291 { 292 PetscErrorCode ierr; 293 PetscBool iascii; 294 PetscViewerFormat format; 295 296 PetscFunctionBegin; 297 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 298 if (iascii) { 299 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 300 if (format == PETSC_VIEWER_ASCII_INFO) { 301 ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 302 } 303 } 304 PetscFunctionReturn(0); 305 } 306 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "MatSolve_SuperLU_Private" 310 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 311 { 312 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 313 const PetscScalar *barray; 314 PetscScalar *xarray; 315 PetscErrorCode ierr; 316 PetscInt info,i,n; 317 PetscReal ferr,berr; 318 static PetscBool cite = PETSC_FALSE; 319 320 PetscFunctionBegin; 321 if (lu->lwork == -1) PetscFunctionReturn(0); 322 ierr = PetscCitationsRegister("@article{superlu99,\n author = {James W. Demmel and Stanley C. Eisenstat and\n John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n title = {A supernodal approach to sparse partial pivoting},\n journal = {SIAM J. Matrix Analysis and Applications},\n year = {1999},\n volume = {20},\n number = {3},\n pages = {720-755}\n}\n",&cite);CHKERRQ(ierr); 323 324 ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 325 lu->B.ncol = 1; /* Set the number of right-hand side */ 326 if (lu->options.Equil && !lu->rhs_dup) { 327 /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 328 ierr = PetscMalloc1(n,&lu->rhs_dup);CHKERRQ(ierr); 329 } 330 if (lu->options.Equil) { 331 /* Copy b into rsh_dup */ 332 ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 333 ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 334 ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 335 barray = lu->rhs_dup; 336 } else { 337 ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 338 } 339 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 340 341 #if defined(PETSC_USE_COMPLEX) 342 #if defined(PETSC_USE_REAL_SINGLE) 343 ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 344 ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 345 #else 346 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 347 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 348 #endif 349 #else 350 ((DNformat*)lu->B.Store)->nzval = (void*)barray; 351 ((DNformat*)lu->X.Store)->nzval = xarray; 352 #endif 353 354 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 355 if (A->factortype == MAT_FACTOR_LU) { 356 #if defined(PETSC_USE_COMPLEX) 357 #if defined(PETSC_USE_REAL_SINGLE) 358 PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 359 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 360 &lu->mem_usage, &lu->stat, &info)); 361 #else 362 PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 363 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 364 &lu->mem_usage, &lu->stat, &info)); 365 #endif 366 #else 367 #if defined(PETSC_USE_REAL_SINGLE) 368 PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 369 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 370 &lu->mem_usage, &lu->stat, &info)); 371 #else 372 PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 373 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 374 &lu->Glu,&lu->mem_usage, &lu->stat, &info)); 375 #endif 376 #endif 377 } else if (A->factortype == MAT_FACTOR_ILU) { 378 #if defined(PETSC_USE_COMPLEX) 379 #if defined(PETSC_USE_REAL_SINGLE) 380 PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 381 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 382 &lu->mem_usage, &lu->stat, &info)); 383 #else 384 PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 385 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 386 &lu->mem_usage, &lu->stat, &info)); 387 #endif 388 #else 389 #if defined(PETSC_USE_REAL_SINGLE) 390 PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 391 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 392 &lu->mem_usage, &lu->stat, &info)); 393 #else 394 PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 395 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 396 &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 397 #endif 398 #endif 399 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 400 if (!lu->options.Equil) { 401 ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 402 } 403 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 404 405 if (!info || info == lu->A.ncol+1) { 406 if (lu->options.IterRefine) { 407 ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 408 ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 409 for (i = 0; i < 1; ++i) { 410 ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 411 } 412 } 413 } else if (info > 0) { 414 if (lu->lwork == -1) { 415 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 416 } else { 417 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 418 } 419 } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 420 421 if (lu->options.PrintStat) { 422 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 423 PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 424 } 425 PetscFunctionReturn(0); 426 } 427 428 #undef __FUNCT__ 429 #define __FUNCT__ "MatSolve_SuperLU" 430 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 431 { 432 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 lu->options.Trans = TRANS; 437 438 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "MatSolveTranspose_SuperLU" 444 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 445 { 446 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 447 PetscErrorCode ierr; 448 449 PetscFunctionBegin; 450 lu->options.Trans = NOTRANS; 451 452 ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 453 PetscFunctionReturn(0); 454 } 455 456 #undef __FUNCT__ 457 #define __FUNCT__ "MatMatSolve_SuperLU" 458 PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 459 { 460 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 461 PetscBool flg; 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 466 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 467 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 468 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 469 lu->options.Trans = TRANS; 470 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 471 PetscFunctionReturn(0); 472 } 473 474 /* 475 Note the r permutation is ignored 476 */ 477 #undef __FUNCT__ 478 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 479 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 480 { 481 Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 482 483 PetscFunctionBegin; 484 lu->flg = DIFFERENT_NONZERO_PATTERN; 485 lu->CleanUpSuperLU = PETSC_TRUE; 486 F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 492 static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 493 { 494 Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 495 496 PetscFunctionBegin; 497 lu->options.ILU_DropTol = dtol; 498 PetscFunctionReturn(0); 499 } 500 501 #undef __FUNCT__ 502 #define __FUNCT__ "MatSuperluSetILUDropTol" 503 /*@ 504 MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 505 Logically Collective on Mat 506 507 Input Parameters: 508 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 509 - dtol - drop tolerance 510 511 Options Database: 512 . -mat_superlu_ilu_droptol <dtol> 513 514 Level: beginner 515 516 References: SuperLU Users' Guide 517 518 .seealso: MatGetFactor() 519 @*/ 520 PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 521 { 522 PetscErrorCode ierr; 523 524 PetscFunctionBegin; 525 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 526 PetscValidLogicalCollectiveReal(F,dtol,2); 527 ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 533 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 534 { 535 PetscFunctionBegin; 536 *type = MATSOLVERSUPERLU; 537 PetscFunctionReturn(0); 538 } 539 540 541 /*MC 542 MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 543 via the external package SuperLU. 544 545 Use ./configure --download-superlu to have PETSc installed with SuperLU 546 547 Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver 548 549 Options Database Keys: 550 + -mat_superlu_equil <FALSE> - Equil (None) 551 . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 552 . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 553 . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 554 . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 555 . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 556 . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 557 . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 558 . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 559 . -mat_superlu_printstat <FALSE> - PrintStat (None) 560 . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 561 . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 562 . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 563 . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 564 . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 565 . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 566 - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 567 568 Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 569 570 Level: beginner 571 572 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 573 M*/ 574 575 #undef __FUNCT__ 576 #define __FUNCT__ "MatGetFactor_seqaij_superlu" 577 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 578 { 579 Mat B; 580 Mat_SuperLU *lu; 581 PetscErrorCode ierr; 582 PetscInt indx,m=A->rmap->n,n=A->cmap->n; 583 PetscBool flg,set; 584 PetscReal real_input; 585 const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 586 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 587 const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 588 589 PetscFunctionBegin; 590 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 591 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 592 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 593 ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 594 595 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 596 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 597 B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 598 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 599 600 B->ops->destroy = MatDestroy_SuperLU; 601 B->ops->view = MatView_SuperLU; 602 B->ops->getdiagonal = MatGetDiagonal_SuperLU; 603 B->factortype = ftype; 604 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 605 B->preallocated = PETSC_TRUE; 606 607 ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 608 609 if (ftype == MAT_FACTOR_LU) { 610 set_default_options(&lu->options); 611 /* Comments from SuperLU_4.0/SRC/dgssvx.c: 612 "Whether or not the system will be equilibrated depends on the 613 scaling of the matrix A, but if equilibration is used, A is 614 overwritten by diag(R)*A*diag(C) and B by diag(R)*B 615 (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 616 We set 'options.Equil = NO' as default because additional space is needed for it. 617 */ 618 lu->options.Equil = NO; 619 } else if (ftype == MAT_FACTOR_ILU) { 620 /* Set the default input options of ilu: */ 621 PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options)); 622 } 623 lu->options.PrintStat = NO; 624 625 /* Initialize the statistics variables. */ 626 PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat)); 627 lu->lwork = 0; /* allocate space internally by system malloc */ 628 629 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 630 ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);CHKERRQ(ierr); 631 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 632 if (flg) lu->options.ColPerm = (colperm_t)indx; 633 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 634 if (flg) lu->options.IterRefine = (IterRefine_t)indx; 635 ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);CHKERRQ(ierr); 636 if (set && flg) lu->options.SymmetricMode = YES; 637 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr); 638 if (flg) lu->options.DiagPivotThresh = (double) real_input; 639 ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);CHKERRQ(ierr); 640 if (set && flg) lu->options.PivotGrowth = YES; 641 ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);CHKERRQ(ierr); 642 if (set && flg) lu->options.ConditionNumber = YES; 643 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 644 if (flg) lu->options.RowPerm = (rowperm_t)indx; 645 ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);CHKERRQ(ierr); 646 if (set && flg) lu->options.ReplaceTinyPivot = YES; 647 ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);CHKERRQ(ierr); 648 if (set && flg) lu->options.PrintStat = YES; 649 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);CHKERRQ(ierr); 650 if (lu->lwork > 0) { 651 /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/ 652 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 653 } else if (lu->lwork != 0 && lu->lwork != -1) { 654 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 655 lu->lwork = 0; 656 } 657 /* ilu options */ 658 ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr); 659 if (flg) lu->options.ILU_DropTol = (double) real_input; 660 ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr); 661 if (flg) lu->options.ILU_FillTol = (double) real_input; 662 ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr); 663 if (flg) lu->options.ILU_FillFactor = (double) real_input; 664 ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);CHKERRQ(ierr); 665 ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 666 if (flg) lu->options.ILU_Norm = (norm_t)indx; 667 ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 668 if (flg) lu->options.ILU_MILU = (milu_t)indx; 669 PetscOptionsEnd(); 670 if (lu->options.Equil == YES) { 671 /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 672 ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 673 } 674 675 /* Allocate spaces (notice sizes are for the transpose) */ 676 ierr = PetscMalloc1(m,&lu->etree);CHKERRQ(ierr); 677 ierr = PetscMalloc1(n,&lu->perm_r);CHKERRQ(ierr); 678 ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr); 679 ierr = PetscMalloc1(n,&lu->R);CHKERRQ(ierr); 680 ierr = PetscMalloc1(m,&lu->C);CHKERRQ(ierr); 681 682 /* create rhs and solution x without allocate space for .Store */ 683 #if defined(PETSC_USE_COMPLEX) 684 #if defined(PETSC_USE_REAL_SINGLE) 685 PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 686 PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 687 #else 688 PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 689 PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 690 #endif 691 #else 692 #if defined(PETSC_USE_REAL_SINGLE) 693 PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 694 PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 695 #else 696 PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 697 PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 698 #endif 699 #endif 700 701 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 702 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 703 B->spptr = lu; 704 705 *F = B; 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "MatSolverPackageRegister_SuperLU" 711 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void) 712 { 713 PetscErrorCode ierr; 714 715 PetscFunctionBegin; 716 ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr); 717 ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720