Actual source code: factor.c

  1: #include <../src/ksp/pc/impls/factor/factor.h>
  2: #include <petsc/private/matimpl.h>

  4: /*
  5:     If an ordering is not yet set and the matrix is available determine a default ordering
  6: */
  7: PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
  8: {
  9:   PetscBool   foundmtype, flg, destroy = PETSC_FALSE;
 10:   const char *prefix;

 12:   PetscFunctionBegin;
 13:   if (pc->pmat) {
 14:     PetscCall(PCGetOptionsPrefix(pc, &prefix));
 15:     PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
 16:     PC_Factor *fact = (PC_Factor *)pc->data;
 17:     PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL));
 18:     if (foundmtype) {
 19:       if (!fact->fact) {
 20:         /* factored matrix is not present at this point, we want to create it during PCSetUp.
 21:            Since this can be called from setfromoptions, we destroy it when we are done with it */
 22:         PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact));
 23:         destroy = PETSC_TRUE;
 24:       }
 25:       if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS);
 26:       if (!fact->fact->assembled) {
 27:         PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg));
 28:         if (!flg) {
 29:           Mat B;
 30:           PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B));
 31:           PetscCall(MatHeaderReplace(fact->fact, &B));
 32:         }
 33:       }
 34:       if (!fact->ordering) {
 35:         PetscBool       canuseordering;
 36:         MatOrderingType otype;

 38:         PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering));
 39:         if (canuseordering) PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype));
 40:         else otype = MATORDERINGEXTERNAL;
 41:         PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering));
 42:       }
 43:       if (destroy) PetscCall(MatDestroy(&fact->fact));
 44:     }
 45:   }
 46:   PetscFunctionReturn(PETSC_SUCCESS);
 47: }

 49: static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag)
 50: {
 51:   PC_Factor *lu = (PC_Factor *)pc->data;

 53:   PetscFunctionBegin;
 54:   lu->reuseordering = flag;
 55:   PetscFunctionReturn(PETSC_SUCCESS);
 56: }

 58: static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag)
 59: {
 60:   PC_Factor *lu = (PC_Factor *)pc->data;

 62:   PetscFunctionBegin;
 63:   lu->reusefill = flag;
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg)
 68: {
 69:   PC_Factor *dir = (PC_Factor *)pc->data;

 71:   PetscFunctionBegin;
 72:   dir->inplace = flg;
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg)
 77: {
 78:   PC_Factor *dir = (PC_Factor *)pc->data;

 80:   PetscFunctionBegin;
 81:   *flg = dir->inplace;
 82:   PetscFunctionReturn(PETSC_SUCCESS);
 83: }

 85: /*@
 86:   PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may
 87:   set the options for that particular factorization object.

 89:   Input Parameter:
 90: . pc - the preconditioner context

 92:   Level: intermediate

 94:   Note:
 95:   After you have called this function (which has to be after the `KSPSetOperators()` or `PCSetOperators()`) you can call `PCFactorGetMatrix()` and then set factor options on that matrix.
 96:   This function raises an error if the requested combination of solver package and matrix type is not supported.

 98:   Developer Note:
 99:   This function should have a name that clearly indicates that this calls `MatGetFactor()` and thus populates the `->factor` field. The `MatSetSolverType` portion of the name
100:   may not add value to the clarity of the purpose of the function and could be removed.

102: .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
103: @*/
104: PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
105: {
106:   PetscFunctionBegin;
108:   PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc));
109:   PetscFunctionReturn(PETSC_SUCCESS);
110: }

112: /*@
113:   PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero

115:   Logically Collective

117:   Input Parameters:
118: + pc   - the preconditioner context
119: - zero - all pivots smaller than this will be considered zero

121:   Options Database Key:
122: . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot

124:   Level: intermediate

126: .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
127: @*/
128: PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero)
129: {
130:   PetscFunctionBegin;
133:   PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero));
134:   PetscFunctionReturn(PETSC_SUCCESS);
135: }

137: /*@
138:   PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
139:   numerical factorization, thus the matrix has nonzero pivots

141:   Logically Collective

143:   Input Parameters:
144: + pc        - the preconditioner context
145: - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS`

147:   Options Database Key:
148: . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types

150:   Level: intermediate

152: .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
153: @*/
154: PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype)
155: {
156:   PetscFunctionBegin;
159:   PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype));
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

163: /*@
164:   PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
165:   numerical factorization, thus the matrix has nonzero pivots

167:   Logically Collective

169:   Input Parameters:
170: + pc          - the preconditioner context
171: - shiftamount - amount of shift or `PETSC_DECIDE` for the default

173:   Options Database Key:
174: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default

176:   Level: intermediate

178: .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
179: @*/
180: PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount)
181: {
182:   PetscFunctionBegin;
185:   PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount));
186:   PetscFunctionReturn(PETSC_SUCCESS);
187: }

189: /*@
190:   PCFactorSetDropTolerance - The preconditioner will use an `PCILU`
191:   based on a drop tolerance.

193:   Logically Collective

195:   Input Parameters:
196: + pc          - the preconditioner context
197: . dt          - the drop tolerance, try from 1.e-10 to .1
198: . dtcol       - tolerance for column pivot, good values [0.1 to 0.01]
199: - maxrowcount - the max number of nonzeros allowed in a row, best value
200:                  depends on the number of nonzeros in row of original matrix

202:   Options Database Key:
203: . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance

205:   Level: intermediate

207:   Note:
208:   There are NO default values for the 3 parameters, you must set them with reasonable values for your
209:   matrix. We don't know how to compute reasonable values.

211: .seealso: [](ch_ksp), `PCILU`
212: @*/
213: PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount)
214: {
215:   PetscFunctionBegin;
219:   PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount));
220:   PetscFunctionReturn(PETSC_SUCCESS);
221: }

223: /*@
224:   PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot

226:   Not Collective

228:   Input Parameter:
229: . pc - the preconditioner context

231:   Output Parameter:
232: . pivot - the tolerance

234:   Level: intermediate

236: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()`
237: @*/
238: PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot)
239: {
240:   PetscFunctionBegin;
242:   PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot));
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }

246: /*@
247:   PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot

249:   Not Collective

251:   Input Parameter:
252: . pc - the preconditioner context

254:   Output Parameter:
255: . shift - how much to shift the diagonal entry

257:   Level: intermediate

259: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
260: @*/
261: PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift)
262: {
263:   PetscFunctionBegin;
265:   PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift));
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: /*@
270:   PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected

272:   Not Collective

274:   Input Parameter:
275: . pc - the preconditioner context

277:   Output Parameter:
278: . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS`

280:   Level: intermediate

282: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
283: @*/
284: PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type)
285: {
286:   PetscFunctionBegin;
288:   PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type));
289:   PetscFunctionReturn(PETSC_SUCCESS);
290: }

292: /*@
293:   PCFactorGetLevels - Gets the number of levels of fill to use.

295:   Logically Collective

297:   Input Parameter:
298: . pc - the preconditioner context

300:   Output Parameter:
301: . levels - number of levels of fill

303:   Level: intermediate

305: .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetLevels()`
306: @*/
307: PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels)
308: {
309:   PetscFunctionBegin;
311:   PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels));
312:   PetscFunctionReturn(PETSC_SUCCESS);
313: }

315: /*@
316:   PCFactorSetLevels - Sets the number of levels of fill to use.

318:   Logically Collective

320:   Input Parameters:
321: + pc     - the preconditioner context
322: - levels - number of levels of fill

324:   Options Database Key:
325: . -pc_factor_levels <levels> - Sets fill level

327:   Level: intermediate

329: .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetLevels()`
330: @*/
331: PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels)
332: {
333:   PetscFunctionBegin;
335:   PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels");
337:   PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels));
338:   PetscFunctionReturn(PETSC_SUCCESS);
339: }

341: /*@
342:   PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
343:   treated as level 0 fill even if there is no non-zero location.

345:   Logically Collective

347:   Input Parameters:
348: + pc  - the preconditioner context
349: - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off

351:   Options Database Key:
352: . -pc_factor_diagonal_fill <bool> - allow the diagonal fill

354:   Note:
355:   Does not apply with 0 fill.

357:   Level: intermediate

359: .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()`
360: @*/
361: PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg)
362: {
363:   PetscFunctionBegin;
365:   PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg));
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: /*@
370:   PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
371:   treated as level 0 fill even if there is no non-zero location.

373:   Logically Collective

375:   Input Parameter:
376: . pc - the preconditioner context

378:   Output Parameter:
379: . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off

381:   Note:
382:   Does not apply with 0 fill.

384:   Level: intermediate

386: .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()`
387: @*/
388: PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg)
389: {
390:   PetscFunctionBegin;
392:   PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg));
393:   PetscFunctionReturn(PETSC_SUCCESS);
394: }

396: /*@
397:   PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal

399:   Logically Collective

401:   Input Parameters:
402: + pc   - the preconditioner context
403: - rtol - diagonal entries smaller than this in absolute value are considered zero

405:   Options Database Key:
406: . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance

408:   Level: intermediate

410: .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftAmount()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
411: @*/
412: PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol)
413: {
414:   PetscFunctionBegin;
417:   PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol));
418:   PetscFunctionReturn(PETSC_SUCCESS);
419: }

421: /*@
422:   PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization

424:   Logically Collective

426:   Input Parameters:
427: + pc    - the preconditioner context
428: - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`

430:   Options Database Key:
431: . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse

433:   Level: intermediate

435:   Note:
436:   The default type is set by searching for available types based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`.
437:   Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver.
438:   For example if one configuration had --download-mumps while a different one had --download-superlu_dist.

440: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverTypeRegister()`,
441:           `MatInitializePackage()`, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `MatSolverTypeGet()`
442: @*/
443: PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype)
444: {
445:   PetscFunctionBegin;
447:   PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype));
448:   PetscFunctionReturn(PETSC_SUCCESS);
449: }

451: /*@
452:   PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization

454:   Not Collective

456:   Input Parameter:
457: . pc - the preconditioner context

459:   Output Parameter:
460: . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`

462:   Level: intermediate

464: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `MATSOLVERSUPERLU`,
465: `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
466: @*/
467: PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype)
468: {
469:   PetscErrorCode (*f)(PC, MatSolverType *);

471:   PetscFunctionBegin;
473:   PetscAssertPointer(stype, 2);
474:   PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f));
475:   if (f) PetscCall((*f)(pc, stype));
476:   else *stype = NULL;
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: /*@
481:   PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
482:   fill = number nonzeros in factor/number nonzeros in original matrix.

484:   Not Collective, each process can expect a different amount of fill

486:   Input Parameters:
487: + pc   - the preconditioner context
488: - fill - amount of expected fill

490:   Options Database Key:
491: . -pc_factor_fill <fill> - Sets fill amount

493:   Level: intermediate

495:   Notes:
496:   For sparse matrix factorizations it is difficult to predict how much
497:   fill to expect. By running with the option -info PETSc will print the
498:   actual amount of fill used; allowing you to set the value accurately for
499:   future runs. Default PETSc uses a value of 5.0

501:   This is ignored for most solver packages

503:   This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels.

505: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()`
506: @*/
507: PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill)
508: {
509:   PetscFunctionBegin;
511:   PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less than 1.0");
512:   PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill));
513:   PetscFunctionReturn(PETSC_SUCCESS);
514: }

516: /*@
517:   PCFactorSetUseInPlace - Tells the preconditioner to do an in-place factorization.

519:   Logically Collective

521:   Input Parameters:
522: + pc  - the preconditioner context
523: - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable

525:   Options Database Key:
526: . -pc_factor_in_place <true,false> - Activate/deactivate in-place factorization

528:   Note:
529:   For dense matrices, this enables the solution of much larger problems.
530:   For sparse matrices the factorization cannot be done truly in-place
531:   so this does not save memory during the factorization, but after the matrix
532:   is factored, the original unfactored matrix is freed, thus recovering that
533:   space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.

535:   `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when
536:   a different matrix is provided for the multiply and the preconditioner in
537:   a call to `KSPSetOperators()`.
538:   This is because the Krylov space methods require an application of the
539:   matrix multiplication, which is not possible here because the matrix has
540:   been factored in-place, replacing the original matrix.

542:   Level: intermediate

544: .seealso: [](ch_ksp), `PC`, `Mat`, `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()`
545: @*/
546: PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg)
547: {
548:   PetscFunctionBegin;
550:   PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg));
551:   PetscFunctionReturn(PETSC_SUCCESS);
552: }

554: /*@
555:   PCFactorGetUseInPlace - Determines if an in-place factorization is being used.

557:   Logically Collective

559:   Input Parameter:
560: . pc - the preconditioner context

562:   Output Parameter:
563: . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable

565:   Level: intermediate

567: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()`
568: @*/
569: PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg)
570: {
571:   PetscFunctionBegin;
573:   PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg));
574:   PetscFunctionReturn(PETSC_SUCCESS);
575: }

577: /*@
578:   PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
579:   be used in the `PCLU`, `PCCHOLESKY`, `PCILU`,  or `PCICC` preconditioners

581:   Logically Collective

583:   Input Parameters:
584: + pc       - the preconditioner context
585: - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM`

587:   Options Database Key:
588: . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine

590:   Level: intermediate

592:   Notes:
593:   Nested dissection is used by default for some of PETSc's sparse matrix formats

595:   For `PCCHOLESKY` and `PCICC` and the `MATSBAIJ` format the only reordering available is natural since only the upper half of the matrix is stored
596:   and reordering this matrix is very expensive.

598:   You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering.

600:   `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others.

602: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM`
603: @*/
604: PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering)
605: {
606:   PetscFunctionBegin;
608:   PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering));
609:   PetscFunctionReturn(PETSC_SUCCESS);
610: }

612: /*@
613:   PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
614:   For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
615:   it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used.

617:   Logically Collective

619:   Input Parameters:
620: + pc    - the preconditioner context
621: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)

623:   Options Database Key:
624: . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance

626:   Level: intermediate

628: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
629: @*/
630: PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol)
631: {
632:   PetscFunctionBegin;
635:   PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol));
636:   PetscFunctionReturn(PETSC_SUCCESS);
637: }

639: /*@
640:   PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
641:   with `MATBAIJ` or `MATSBAIJ` matrices

643:   Logically Collective

645:   Input Parameters:
646: + pc    - the preconditioner context
647: - pivot - `PETSC_TRUE` or `PETSC_FALSE`

649:   Options Database Key:
650: . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ`

652:   Level: intermediate

654: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
655: @*/
656: PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot)
657: {
658:   PetscFunctionBegin;
661:   PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot));
662:   PetscFunctionReturn(PETSC_SUCCESS);
663: }

665: /*@
666:   PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
667:   this causes later ones to use the fill ratio computed in the initial factorization.

669:   Logically Collective

671:   Input Parameters:
672: + pc   - the preconditioner context
673: - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`

675:   Options Database Key:
676: . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`

678:   Level: intermediate

680: .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()`
681: @*/
682: PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag)
683: {
684:   PetscFunctionBegin;
687:   PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag));
688:   PetscFunctionReturn(PETSC_SUCCESS);
689: }

691: PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype)
692: {
693:   PC_Factor *fact = (PC_Factor *)pc->data;

695:   PetscFunctionBegin;
696:   PetscCall(MatFactorInfoInitialize(&fact->info));
697:   fact->factortype           = ftype;
698:   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
699:   fact->info.shiftamount     = 100.0 * PETSC_MACHINE_EPSILON;
700:   fact->info.zeropivot       = 100.0 * PETSC_MACHINE_EPSILON;
701:   fact->info.pivotinblocks   = 1.0;
702:   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;

704:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor));
705:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor));
706:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor));
707:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor));
708:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor));
709:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor));
710:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor));
711:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor));
712:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor));
713:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor));
714:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor));
715:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor));
716:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor));
717:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor));
718:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor));
719:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor));
720:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor));
721:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor));
722:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor));
723:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor));
724:   PetscFunctionReturn(PETSC_SUCCESS);
725: }

727: PetscErrorCode PCFactorClearComposedFunctions(PC pc)
728: {
729:   PetscFunctionBegin;
730:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL));
731:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL));
732:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL));
733:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL));
734:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL));
735:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL));
736:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL));
737:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL));
738:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL));
739:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL));
740:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL));
741:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL));
742:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL));
743:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL));
744:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL));
745:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL));
746:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL));
747:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL));
748:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL));
749:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL));
750:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL));
751:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL));
752:   PetscFunctionReturn(PETSC_SUCCESS);
753: }