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