1 2 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 3 4 PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC pc) 5 { 6 PC_Factor *icc = (PC_Factor *)pc->data; 7 8 PetscFunctionBegin; 9 PetscCheck(pc->pmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "You can only call this routine after the matrix object has been provided to the solver, for example with KSPSetOperators() or SNESSetJacobian()"); 10 if (!((PC_Factor *)icc)->fact) PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)icc)->solvertype, ((PC_Factor *)icc)->factortype, &((PC_Factor *)icc)->fact)); 11 PetscCheck(((PC_Factor *)icc)->fact, PetscObjectComm((PetscObject)pc->pmat), PETSC_ERR_SUP, "MatFactor type %s not supported by this matrix instance of type %s and solver type %s. Rerun with \"-info :mat | grep MatGetFactor_\" for additional information.", 12 MatFactorTypes[((PC_Factor *)icc)->factortype], ((PetscObject)pc->pmat)->type_name, ((PC_Factor *)icc)->solvertype); 13 PetscFunctionReturn(PETSC_SUCCESS); 14 } 15 16 PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc, PetscReal z) 17 { 18 PC_Factor *ilu = (PC_Factor *)pc->data; 19 20 PetscFunctionBegin; 21 ilu->info.zeropivot = z; 22 PetscFunctionReturn(PETSC_SUCCESS); 23 } 24 25 PetscErrorCode PCFactorSetShiftType_Factor(PC pc, MatFactorShiftType shifttype) 26 { 27 PC_Factor *dir = (PC_Factor *)pc->data; 28 29 PetscFunctionBegin; 30 if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 31 else { 32 dir->info.shifttype = (PetscReal)shifttype; 33 if ((shifttype == MAT_SHIFT_NONZERO || shifttype == MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) { dir->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ } 34 } 35 PetscFunctionReturn(PETSC_SUCCESS); 36 } 37 38 PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc, PetscReal shiftamount) 39 { 40 PC_Factor *dir = (PC_Factor *)pc->data; 41 42 PetscFunctionBegin; 43 if (shiftamount == (PetscReal)PETSC_DECIDE) dir->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 44 else dir->info.shiftamount = shiftamount; 45 PetscFunctionReturn(PETSC_SUCCESS); 46 } 47 48 PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount) 49 { 50 PC_Factor *ilu = (PC_Factor *)pc->data; 51 52 PetscFunctionBegin; 53 if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor *)ilu)->info.dt != dt || ((PC_Factor *)ilu)->info.dtcol != dtcol || ((PC_Factor *)ilu)->info.dtcount != dtcount)) { 54 SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change tolerance after use"); 55 } 56 ilu->info.usedt = PETSC_TRUE; 57 ilu->info.dt = dt; 58 ilu->info.dtcol = dtcol; 59 ilu->info.dtcount = dtcount; 60 ilu->info.fill = PETSC_DEFAULT; 61 PetscFunctionReturn(PETSC_SUCCESS); 62 } 63 64 PetscErrorCode PCFactorSetFill_Factor(PC pc, PetscReal fill) 65 { 66 PC_Factor *dir = (PC_Factor *)pc->data; 67 68 PetscFunctionBegin; 69 dir->info.fill = fill; 70 PetscFunctionReturn(PETSC_SUCCESS); 71 } 72 73 PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc, MatOrderingType ordering) 74 { 75 PC_Factor *dir = (PC_Factor *)pc->data; 76 PetscBool flg; 77 78 PetscFunctionBegin; 79 if (!pc->setupcalled) { 80 PetscCall(PetscFree(dir->ordering)); 81 PetscCall(PetscStrallocpy(ordering, (char **)&dir->ordering)); 82 } else { 83 PetscCall(PetscStrcmp(dir->ordering, ordering, &flg)); 84 PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change ordering after use"); 85 } 86 PetscFunctionReturn(PETSC_SUCCESS); 87 } 88 89 PetscErrorCode PCFactorGetLevels_Factor(PC pc, PetscInt *levels) 90 { 91 PC_Factor *ilu = (PC_Factor *)pc->data; 92 93 PetscFunctionBegin; 94 *levels = ilu->info.levels; 95 PetscFunctionReturn(PETSC_SUCCESS); 96 } 97 98 PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc, PetscReal *pivot) 99 { 100 PC_Factor *ilu = (PC_Factor *)pc->data; 101 102 PetscFunctionBegin; 103 *pivot = ilu->info.zeropivot; 104 PetscFunctionReturn(PETSC_SUCCESS); 105 } 106 107 PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc, PetscReal *shift) 108 { 109 PC_Factor *ilu = (PC_Factor *)pc->data; 110 111 PetscFunctionBegin; 112 *shift = ilu->info.shiftamount; 113 PetscFunctionReturn(PETSC_SUCCESS); 114 } 115 116 PetscErrorCode PCFactorGetShiftType_Factor(PC pc, MatFactorShiftType *type) 117 { 118 PC_Factor *ilu = (PC_Factor *)pc->data; 119 120 PetscFunctionBegin; 121 *type = (MatFactorShiftType)(int)ilu->info.shifttype; 122 PetscFunctionReturn(PETSC_SUCCESS); 123 } 124 125 PetscErrorCode PCFactorSetLevels_Factor(PC pc, PetscInt levels) 126 { 127 PC_Factor *ilu = (PC_Factor *)pc->data; 128 129 PetscFunctionBegin; 130 if (!pc->setupcalled) ilu->info.levels = levels; 131 else if (ilu->info.levels != levels) { 132 PetscUseTypeMethod(pc, reset); /* remove previous factored matrices */ 133 pc->setupcalled = 0; /* force a complete rebuild of preconditioner factored matrices */ 134 ilu->info.levels = levels; 135 } else PetscCheck(!ilu->info.usedt, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change levels after use with ILUdt"); 136 PetscFunctionReturn(PETSC_SUCCESS); 137 } 138 139 PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc, PetscBool flg) 140 { 141 PC_Factor *dir = (PC_Factor *)pc->data; 142 143 PetscFunctionBegin; 144 dir->info.diagonal_fill = (PetscReal)flg; 145 PetscFunctionReturn(PETSC_SUCCESS); 146 } 147 148 PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc, PetscBool *flg) 149 { 150 PC_Factor *dir = (PC_Factor *)pc->data; 151 152 PetscFunctionBegin; 153 *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE; 154 PetscFunctionReturn(PETSC_SUCCESS); 155 } 156 157 PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc, PetscBool pivot) 158 { 159 PC_Factor *dir = (PC_Factor *)pc->data; 160 161 PetscFunctionBegin; 162 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 163 PetscFunctionReturn(PETSC_SUCCESS); 164 } 165 166 PetscErrorCode PCFactorGetMatrix_Factor(PC pc, Mat *mat) 167 { 168 PC_Factor *ilu = (PC_Factor *)pc->data; 169 170 PetscFunctionBegin; 171 PetscCheck(ilu->fact, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 172 *mat = ilu->fact; 173 PetscFunctionReturn(PETSC_SUCCESS); 174 } 175 176 /* allow access to preallocation information */ 177 #include <petsc/private/matimpl.h> 178 179 PetscErrorCode PCFactorSetMatSolverType_Factor(PC pc, MatSolverType stype) 180 { 181 PC_Factor *lu = (PC_Factor *)pc->data; 182 183 PetscFunctionBegin; 184 if (lu->fact && lu->fact->assembled) { 185 MatSolverType ltype; 186 PetscBool flg; 187 PetscCall(MatFactorGetSolverType(lu->fact, <ype)); 188 PetscCall(PetscStrcmp(stype, ltype, &flg)); 189 PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change solver matrix package from %s to %s after PC has been setup or used", ltype, stype); 190 } 191 192 PetscCall(PetscFree(lu->solvertype)); 193 PetscCall(PetscStrallocpy(stype, &lu->solvertype)); 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 PetscErrorCode PCFactorGetMatSolverType_Factor(PC pc, MatSolverType *stype) 198 { 199 PC_Factor *lu = (PC_Factor *)pc->data; 200 201 PetscFunctionBegin; 202 *stype = lu->solvertype; 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc, PetscReal dtcol) 207 { 208 PC_Factor *dir = (PC_Factor *)pc->data; 209 210 PetscFunctionBegin; 211 PetscCheck(dtcol >= 0.0 && dtcol <= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Column pivot tolerance is %g must be between 0 and 1", (double)dtcol); 212 dir->info.dtcol = dtcol; 213 PetscFunctionReturn(PETSC_SUCCESS); 214 } 215 216 PetscErrorCode PCSetFromOptions_Factor(PC pc, PetscOptionItems *PetscOptionsObject) 217 { 218 PC_Factor *factor = (PC_Factor *)pc->data; 219 PetscBool flg, set; 220 char tname[256], solvertype[64]; 221 PetscFunctionList ordlist; 222 PetscEnum etmp; 223 PetscBool inplace; 224 225 PetscFunctionBegin; 226 PetscCall(PCFactorGetUseInPlace(pc, &inplace)); 227 PetscCall(PetscOptionsBool("-pc_factor_in_place", "Form factored matrix in the same memory as the matrix", "PCFactorSetUseInPlace", inplace, &flg, &set)); 228 if (set) PetscCall(PCFactorSetUseInPlace(pc, flg)); 229 PetscCall(PetscOptionsReal("-pc_factor_fill", "Expected non-zeros in factored matrix", "PCFactorSetFill", ((PC_Factor *)factor)->info.fill, &((PC_Factor *)factor)->info.fill, NULL)); 230 231 PetscCall(PetscOptionsEnum("-pc_factor_shift_type", "Type of shift to add to diagonal", "PCFactorSetShiftType", MatFactorShiftTypes, (PetscEnum)(int)((PC_Factor *)factor)->info.shifttype, &etmp, &flg)); 232 if (flg) PetscCall(PCFactorSetShiftType(pc, (MatFactorShiftType)etmp)); 233 PetscCall(PetscOptionsReal("-pc_factor_shift_amount", "Shift added to diagonal", "PCFactorSetShiftAmount", ((PC_Factor *)factor)->info.shiftamount, &((PC_Factor *)factor)->info.shiftamount, NULL)); 234 235 PetscCall(PetscOptionsReal("-pc_factor_zeropivot", "Pivot is considered zero if less than", "PCFactorSetZeroPivot", ((PC_Factor *)factor)->info.zeropivot, &((PC_Factor *)factor)->info.zeropivot, NULL)); 236 PetscCall(PetscOptionsReal("-pc_factor_column_pivot", "Column pivot tolerance (used only for some factorization)", "PCFactorSetColumnPivot", ((PC_Factor *)factor)->info.dtcol, &((PC_Factor *)factor)->info.dtcol, &flg)); 237 238 PetscCall(PetscOptionsBool("-pc_factor_pivot_in_blocks", "Pivot inside matrix dense blocks for BAIJ and SBAIJ", "PCFactorSetPivotInBlocks", ((PC_Factor *)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE, &flg, &set)); 239 if (set) PetscCall(PCFactorSetPivotInBlocks(pc, flg)); 240 241 PetscCall(PetscOptionsBool("-pc_factor_reuse_fill", "Use fill from previous factorization", "PCFactorSetReuseFill", PETSC_FALSE, &flg, &set)); 242 if (set) PetscCall(PCFactorSetReuseFill(pc, flg)); 243 PetscCall(PetscOptionsBool("-pc_factor_reuse_ordering", "Reuse ordering from previous factorization", "PCFactorSetReuseOrdering", PETSC_FALSE, &flg, &set)); 244 if (set) PetscCall(PCFactorSetReuseOrdering(pc, flg)); 245 246 PetscCall(PetscOptionsDeprecated("-pc_factor_mat_solver_package", "-pc_factor_mat_solver_type", "3.9", NULL)); 247 PetscCall(PetscOptionsString("-pc_factor_mat_solver_type", "Specific direct solver to use", "MatGetFactor", ((PC_Factor *)factor)->solvertype, solvertype, sizeof(solvertype), &flg)); 248 if (flg) PetscCall(PCFactorSetMatSolverType(pc, solvertype)); 249 PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 250 PetscCall(MatGetOrderingList(&ordlist)); 251 PetscCall(PetscOptionsFList("-pc_factor_mat_ordering_type", "Reordering to reduce nonzeros in factored matrix", "PCFactorSetMatOrderingType", ordlist, ((PC_Factor *)factor)->ordering, tname, sizeof(tname), &flg)); 252 if (flg) PetscCall(PCFactorSetMatOrderingType(pc, tname)); 253 PetscFunctionReturn(PETSC_SUCCESS); 254 } 255 256 PetscErrorCode PCView_Factor(PC pc, PetscViewer viewer) 257 { 258 PC_Factor *factor = (PC_Factor *)pc->data; 259 PetscBool isstring, iascii, canuseordering; 260 MatInfo info; 261 MatOrderingType ordering; 262 PetscViewerFormat format; 263 264 PetscFunctionBegin; 265 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 266 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 267 if (iascii) { 268 if (factor->inplace) { 269 PetscCall(PetscViewerASCIIPrintf(viewer, " in-place factorization\n")); 270 } else { 271 PetscCall(PetscViewerASCIIPrintf(viewer, " out-of-place factorization\n")); 272 } 273 274 if (factor->reusefill) PetscCall(PetscViewerASCIIPrintf(viewer, " Reusing fill from past factorization\n")); 275 if (factor->reuseordering) PetscCall(PetscViewerASCIIPrintf(viewer, " Reusing reordering from past factorization\n")); 276 if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) { 277 if (factor->info.dt > 0) { 278 PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance %g\n", (double)factor->info.dt)); 279 PetscCall(PetscViewerASCIIPrintf(viewer, " max nonzeros per row %" PetscInt_FMT "\n", (PetscInt)factor->info.dtcount)); 280 PetscCall(PetscViewerASCIIPrintf(viewer, " column permutation tolerance %g\n", (double)factor->info.dtcol)); 281 } else if (factor->info.levels == 1) { 282 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " level of fill\n", (PetscInt)factor->info.levels)); 283 } else { 284 PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " levels of fill\n", (PetscInt)factor->info.levels)); 285 } 286 } 287 288 PetscCall(PetscViewerASCIIPrintf(viewer, " tolerance for zero pivot %g\n", (double)factor->info.zeropivot)); 289 if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */ 290 PetscCall(PetscViewerASCIIPrintf(viewer, " using %s [%s]\n", MatFactorShiftTypesDetail[(int)factor->info.shifttype], MatFactorShiftTypes[(int)factor->info.shifttype])); 291 } 292 293 if (factor->fact) { 294 PetscCall(MatFactorGetCanUseOrdering(factor->fact, &canuseordering)); 295 if (!canuseordering) ordering = MATORDERINGEXTERNAL; 296 else ordering = factor->ordering; 297 PetscCall(PetscViewerASCIIPrintf(viewer, " matrix ordering: %s\n", ordering)); 298 if (!factor->fact->assembled) { 299 PetscCall(PetscViewerASCIIPrintf(viewer, " matrix solver type: %s\n", factor->fact->solvertype)); 300 PetscCall(PetscViewerASCIIPrintf(viewer, " matrix not yet factored; no additional information available\n")); 301 } else { 302 PetscCall(MatGetInfo(factor->fact, MAT_LOCAL, &info)); 303 PetscCall(PetscViewerASCIIPrintf(viewer, " factor fill ratio given %g, needed %g\n", (double)info.fill_ratio_given, (double)info.fill_ratio_needed)); 304 PetscCall(PetscViewerASCIIPrintf(viewer, " Factored matrix follows:\n")); 305 PetscCall(PetscViewerASCIIPushTab(viewer)); 306 PetscCall(PetscViewerASCIIPushTab(viewer)); 307 PetscCall(PetscViewerASCIIPushTab(viewer)); 308 PetscCall(PetscViewerGetFormat(viewer, &format)); 309 PetscCall(PetscViewerPushFormat(viewer, format != PETSC_VIEWER_ASCII_INFO_DETAIL ? PETSC_VIEWER_ASCII_INFO : PETSC_VIEWER_ASCII_INFO_DETAIL)); 310 PetscCall(MatView(factor->fact, viewer)); 311 PetscCall(PetscViewerPopFormat(viewer)); 312 PetscCall(PetscViewerASCIIPopTab(viewer)); 313 PetscCall(PetscViewerASCIIPopTab(viewer)); 314 PetscCall(PetscViewerASCIIPopTab(viewer)); 315 } 316 } 317 318 } else if (isstring) { 319 MatFactorType t; 320 PetscCall(MatGetFactorType(factor->fact, &t)); 321 if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) PetscCall(PetscViewerStringSPrintf(viewer, " lvls=%" PetscInt_FMT ",order=%s", (PetscInt)factor->info.levels, factor->ordering)); 322 } 323 PetscFunctionReturn(PETSC_SUCCESS); 324 } 325