1 2 /* 3 Defines a direct factorization preconditioner for any Mat implementation 4 Note: this need not be consided a preconditioner since it supplies 5 a direct solver. 6 */ 7 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 8 9 typedef struct { 10 PC_Factor hdr; 11 IS row, col; /* index sets used for reordering */ 12 } PC_Cholesky; 13 14 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems *PetscOptionsObject) 15 { 16 PetscFunctionBegin; 17 PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options"); 18 PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject)); 19 PetscOptionsHeadEnd(); 20 PetscFunctionReturn(PETSC_SUCCESS); 21 } 22 23 static PetscErrorCode PCSetUp_Cholesky(PC pc) 24 { 25 PetscBool flg; 26 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 27 MatSolverType stype; 28 MatFactorError err; 29 const char *prefix; 30 31 PetscFunctionBegin; 32 pc->failedreason = PC_NOERROR; 33 if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill; 34 35 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 36 PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 37 38 PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 39 if (dir->hdr.inplace) { 40 MatFactorType ftype; 41 42 PetscCall(MatGetFactorType(pc->pmat, &ftype)); 43 if (ftype == MAT_FACTOR_NONE) { 44 if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row)); 45 PetscCall(ISDestroy(&dir->col)); 46 /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 47 PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 48 PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 49 if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 50 PetscCall(ISDestroy(&dir->col)); 51 } 52 PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 53 PetscCall(MatFactorGetError(pc->pmat, &err)); 54 if (err) { /* Factor() fails */ 55 pc->failedreason = (PCFailedReason)err; 56 PetscFunctionReturn(PETSC_SUCCESS); 57 } 58 } 59 ((PC_Factor *)dir)->fact = pc->pmat; 60 } else { 61 MatInfo info; 62 63 if (!pc->setupcalled) { 64 PetscBool canuseordering; 65 66 PetscCall(PCFactorSetUpMatSolverType(pc)); 67 PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 68 if (canuseordering) { 69 PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 70 PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 71 /* check if dir->row == dir->col */ 72 if (dir->row) { 73 PetscCall(ISEqual(dir->row, dir->col, &flg)); 74 PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal"); 75 } 76 PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */ 77 78 flg = PETSC_FALSE; 79 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 80 if (flg) { 81 PetscReal tol = 1.e-10; 82 PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 83 PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 84 } 85 } 86 PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 87 PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 88 dir->hdr.actualfill = info.fill_ratio_needed; 89 } else if (pc->flag != SAME_NONZERO_PATTERN) { 90 if (!dir->hdr.reuseordering) { 91 PetscBool canuseordering; 92 93 PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 94 PetscCall(PCFactorSetUpMatSolverType(pc)); 95 PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 96 if (canuseordering) { 97 PetscCall(ISDestroy(&dir->row)); 98 PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 99 PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 100 PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */ 101 102 flg = PETSC_FALSE; 103 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 104 if (flg) { 105 PetscReal tol = 1.e-10; 106 PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 107 PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 108 } 109 } 110 } 111 PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 112 PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 113 dir->hdr.actualfill = info.fill_ratio_needed; 114 } else { 115 PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 116 if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 117 PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 118 pc->failedreason = PC_NOERROR; 119 } 120 } 121 PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 122 if (err) { /* FactorSymbolic() fails */ 123 pc->failedreason = (PCFailedReason)err; 124 PetscFunctionReturn(PETSC_SUCCESS); 125 } 126 127 PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 128 PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 129 if (err) { /* FactorNumeric() fails */ 130 pc->failedreason = (PCFailedReason)err; 131 } 132 } 133 134 PetscCall(PCFactorGetMatSolverType(pc, &stype)); 135 if (!stype) { 136 MatSolverType solverpackage; 137 PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 138 PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 139 } 140 PetscFunctionReturn(PETSC_SUCCESS); 141 } 142 143 static PetscErrorCode PCReset_Cholesky(PC pc) 144 { 145 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 146 147 PetscFunctionBegin; 148 if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 149 PetscCall(ISDestroy(&dir->row)); 150 PetscCall(ISDestroy(&dir->col)); 151 PetscFunctionReturn(PETSC_SUCCESS); 152 } 153 154 static PetscErrorCode PCDestroy_Cholesky(PC pc) 155 { 156 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 157 158 PetscFunctionBegin; 159 PetscCall(PCReset_Cholesky(pc)); 160 PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 161 PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 162 PetscCall(PCFactorClearComposedFunctions(pc)); 163 PetscCall(PetscFree(pc->data)); 164 PetscFunctionReturn(PETSC_SUCCESS); 165 } 166 167 static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) 168 { 169 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 170 171 PetscFunctionBegin; 172 if (dir->hdr.inplace) { 173 PetscCall(MatSolve(pc->pmat, x, y)); 174 } else { 175 PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 176 } 177 PetscFunctionReturn(PETSC_SUCCESS); 178 } 179 180 static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) 181 { 182 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 183 184 PetscFunctionBegin; 185 if (dir->hdr.inplace) { 186 PetscCall(MatMatSolve(pc->pmat, X, Y)); 187 } else { 188 PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 189 } 190 PetscFunctionReturn(PETSC_SUCCESS); 191 } 192 193 static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) 194 { 195 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 196 197 PetscFunctionBegin; 198 if (dir->hdr.inplace) { 199 PetscCall(MatForwardSolve(pc->pmat, x, y)); 200 } else { 201 PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y)); 202 } 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) 207 { 208 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 209 210 PetscFunctionBegin; 211 if (dir->hdr.inplace) { 212 PetscCall(MatBackwardSolve(pc->pmat, x, y)); 213 } else { 214 PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y)); 215 } 216 PetscFunctionReturn(PETSC_SUCCESS); 217 } 218 219 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) 220 { 221 PC_Cholesky *dir = (PC_Cholesky *)pc->data; 222 223 PetscFunctionBegin; 224 if (dir->hdr.inplace) { 225 PetscCall(MatSolveTranspose(pc->pmat, x, y)); 226 } else { 227 PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 228 } 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 232 /*@ 233 PCFactorSetReuseOrdering - When similar matrices are factored, this 234 causes the ordering computed in the first factor to be used for all 235 following factors. 236 237 Logically Collective 238 239 Input Parameters: 240 + pc - the preconditioner context 241 - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 242 243 Options Database Key: 244 . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 245 246 Level: intermediate 247 248 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()` 249 @*/ 250 PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) 251 { 252 PetscFunctionBegin; 253 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 254 PetscValidLogicalCollectiveBool(pc, flag, 2); 255 PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag)); 256 PetscFunctionReturn(PETSC_SUCCESS); 257 } 258 259 /*MC 260 PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 261 262 Options Database Keys: 263 + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 264 . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 265 . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 266 . -pc_factor_fill <fill> - Sets fill amount 267 . -pc_factor_in_place - Activates in-place factorization 268 - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 269 270 Level: beginner 271 272 Notes: 273 Not all options work for all matrix formats 274 275 Usually this will compute an "exact" solution in one iteration and does 276 not need a Krylov method (i.e. you can use -ksp_type preonly, or 277 `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 278 279 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 280 `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 281 `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 282 `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()` 283 M*/ 284 285 PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 286 { 287 PC_Cholesky *dir; 288 289 PetscFunctionBegin; 290 PetscCall(PetscNew(&dir)); 291 pc->data = (void *)dir; 292 PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY)); 293 294 ((PC_Factor *)dir)->info.fill = 5.0; 295 296 pc->ops->destroy = PCDestroy_Cholesky; 297 pc->ops->reset = PCReset_Cholesky; 298 pc->ops->apply = PCApply_Cholesky; 299 pc->ops->matapply = PCMatApply_Cholesky; 300 pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 301 pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 302 pc->ops->applytranspose = PCApplyTranspose_Cholesky; 303 pc->ops->setup = PCSetUp_Cholesky; 304 pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 305 pc->ops->view = PCView_Factor; 306 pc->ops->applyrichardson = NULL; 307 PetscFunctionReturn(PETSC_SUCCESS); 308 } 309