Actual source code: precon.c

  1: /*
  2:     The PC (preconditioner) interface routines, callable by users.
  3: */
  4: #include <petsc/private/pcimpl.h>
  5: #include <petscdm.h>

  7: /* Logging support */
  8: PetscClassId  PC_CLASSID;
  9: PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
 10: PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
 11: PetscInt      PetscMGLevelId;
 12: PetscLogStage PCMPIStage;

 14: PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
 15: {
 16:   PetscMPIInt size;
 17:   PetscBool   hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;

 19:   PetscFunctionBegin;
 20:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
 21:   if (pc->pmat) {
 22:     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
 23:     PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
 24:     if (size == 1) {
 25:       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
 26:       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
 27:       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
 28:       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
 29:       if (flg1 && (!flg2 || (set && flg3))) {
 30:         *type = PCICC;
 31:       } else if (flg2) {
 32:         *type = PCILU;
 33:       } else if (isnormal) {
 34:         *type = PCNONE;
 35:       } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
 36:         *type = PCBJACOBI;
 37:       } else if (hasopsolve) {
 38:         *type = PCMAT;
 39:       } else {
 40:         *type = PCNONE;
 41:       }
 42:     } else {
 43:       if (hasopblock) {
 44:         *type = PCBJACOBI;
 45:       } else if (hasopsolve) {
 46:         *type = PCMAT;
 47:       } else {
 48:         *type = PCNONE;
 49:       }
 50:     }
 51:   } else *type = NULL;
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }

 55: /*@
 56:   PCReset - Resets a `PC` context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s

 58:   Collective

 60:   Input Parameter:
 61: . pc - the preconditioner context

 63:   Level: developer

 65:   Note:
 66:   This allows a `PC` to be reused for a different sized linear system but using the same options that have been previously set in `pc`

 68: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
 69: @*/
 70: PetscErrorCode PCReset(PC pc)
 71: {
 72:   PetscFunctionBegin;
 74:   PetscTryTypeMethod(pc, reset);
 75:   PetscCall(VecDestroy(&pc->diagonalscaleright));
 76:   PetscCall(VecDestroy(&pc->diagonalscaleleft));
 77:   PetscCall(MatDestroy(&pc->pmat));
 78:   PetscCall(MatDestroy(&pc->mat));

 80:   pc->setupcalled = 0;
 81:   PetscFunctionReturn(PETSC_SUCCESS);
 82: }

 84: /*@
 85:   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.

 87:   Collective

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

 92:   Level: developer

 94: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
 95: @*/
 96: PetscErrorCode PCDestroy(PC *pc)
 97: {
 98:   PetscFunctionBegin;
 99:   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
101:   if (--((PetscObject)*pc)->refct > 0) {
102:     *pc = NULL;
103:     PetscFunctionReturn(PETSC_SUCCESS);
104:   }

106:   PetscCall(PCReset(*pc));

108:   /* if memory was published with SAWs then destroy it */
109:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
110:   PetscTryTypeMethod(*pc, destroy);
111:   PetscCall(DMDestroy(&(*pc)->dm));
112:   PetscCall(PetscHeaderDestroy(pc));
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: }

116: /*@
117:   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
118:   scaling as needed by certain time-stepping codes.

120:   Logically Collective

122:   Input Parameter:
123: . pc - the preconditioner context

125:   Output Parameter:
126: . flag - `PETSC_TRUE` if it applies the scaling

128:   Level: developer

130:   Note:
131:   If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning,

133:   $$
134:   \begin{align*}
135:   D M A D^{-1} y = D M b  \\
136:   D A M D^{-1} z = D b.
137:   \end{align*}
138:   $$

140: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
141: @*/
142: PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
143: {
144:   PetscFunctionBegin;
146:   PetscAssertPointer(flag, 2);
147:   *flag = pc->diagonalscale;
148:   PetscFunctionReturn(PETSC_SUCCESS);
149: }

151: /*@
152:   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
153:   scaling as needed by certain time-stepping codes.

155:   Logically Collective

157:   Input Parameters:
158: + pc - the preconditioner context
159: - s  - scaling vector

161:   Level: intermediate

163:   Notes:
164:   The system solved via the Krylov method is, for left and right preconditioning,
165:   $$
166:   \begin{align*}
167:   D M A D^{-1} y = D M b \\
168:   D A M D^{-1} z = D b.
169:   \end{align*}
170:   $$

172:   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.

174: .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
175: @*/
176: PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
177: {
178:   PetscFunctionBegin;
181:   pc->diagonalscale = PETSC_TRUE;

183:   PetscCall(PetscObjectReference((PetscObject)s));
184:   PetscCall(VecDestroy(&pc->diagonalscaleleft));

186:   pc->diagonalscaleleft = s;

188:   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
189:   PetscCall(VecCopy(s, pc->diagonalscaleright));
190:   PetscCall(VecReciprocal(pc->diagonalscaleright));
191:   PetscFunctionReturn(PETSC_SUCCESS);
192: }

194: /*@
195:   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.

197:   Logically Collective

199:   Input Parameters:
200: + pc  - the preconditioner context
201: . in  - input vector
202: - out - scaled vector (maybe the same as in)

204:   Level: intermediate

206:   Notes:
207:   The system solved via the Krylov method is, for left and right preconditioning,

209:   $$
210:   \begin{align*}
211:   D M A D^{-1} y = D M b  \\
212:   D A M D^{-1} z = D b.
213:   \end{align*}
214:   $$

216:   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.

218:   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`

220: .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `MatDiagonalScale()`
221: @*/
222: PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
223: {
224:   PetscFunctionBegin;
228:   if (pc->diagonalscale) {
229:     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
230:   } else if (in != out) {
231:     PetscCall(VecCopy(in, out));
232:   }
233:   PetscFunctionReturn(PETSC_SUCCESS);
234: }

236: /*@
237:   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.

239:   Logically Collective

241:   Input Parameters:
242: + pc  - the preconditioner context
243: . in  - input vector
244: - out - scaled vector (maybe the same as in)

246:   Level: intermediate

248:   Notes:
249:   The system solved via the Krylov method is, for left and right preconditioning,

251:   $$
252:   \begin{align*}
253:   D M A D^{-1} y = D M b  \\
254:   D A M D^{-1} z = D b.
255:   \end{align*}
256:   $$

258:   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.

260:   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`

262: .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `MatDiagonalScale()`
263: @*/
264: PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
265: {
266:   PetscFunctionBegin;
270:   if (pc->diagonalscale) {
271:     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
272:   } else if (in != out) {
273:     PetscCall(VecCopy(in, out));
274:   }
275:   PetscFunctionReturn(PETSC_SUCCESS);
276: }

278: /*@
279:   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
280:   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
281:   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.

283:   Logically Collective

285:   Input Parameters:
286: + pc  - the preconditioner context
287: - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)

289:   Options Database Key:
290: . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator

292:   Level: intermediate

294:   Note:
295:   For the common case in which the linear system matrix and the matrix used to construct the
296:   preconditioner are identical, this routine has no affect.

298: .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`,
299:           `KSPSetOperators()`, `PCSetOperators()`
300: @*/
301: PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
302: {
303:   PetscFunctionBegin;
305:   pc->useAmat = flg;
306:   PetscFunctionReturn(PETSC_SUCCESS);
307: }

309: /*@
310:   PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected.

312:   Logically Collective

314:   Input Parameters:
315: + pc  - iterative context obtained from `PCCreate()`
316: - flg - `PETSC_TRUE` indicates you want the error generated

318:   Level: advanced

320:   Notes:
321:   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
322:   to determine if it has converged or failed. Or use -ksp_error_if_not_converged to cause the program to terminate as soon as lack of convergence is
323:   detected.

325:   This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s

327: .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
328: @*/
329: PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
330: {
331:   PetscFunctionBegin;
334:   pc->erroriffailure = flg;
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

338: /*@
339:   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
340:   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
341:   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.

343:   Logically Collective

345:   Input Parameter:
346: . pc - the preconditioner context

348:   Output Parameter:
349: . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)

351:   Level: intermediate

353:   Note:
354:   For the common case in which the linear system matrix and the matrix used to construct the
355:   preconditioner are identical, this routine is does nothing.

357: .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
358: @*/
359: PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
360: {
361:   PetscFunctionBegin;
363:   *flg = pc->useAmat;
364:   PetscFunctionReturn(PETSC_SUCCESS);
365: }

367: /*@
368:   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has

370:   Collective

372:   Input Parameters:
373: + pc    - the `PC`
374: - level - the nest level

376:   Level: developer

378: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
379: @*/
380: PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
381: {
382:   PetscFunctionBegin;
385:   pc->kspnestlevel = level;
386:   PetscFunctionReturn(PETSC_SUCCESS);
387: }

389: /*@
390:   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has

392:   Not Collective

394:   Input Parameter:
395: . pc - the `PC`

397:   Output Parameter:
398: . level - the nest level

400:   Level: developer

402: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
403: @*/
404: PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
405: {
406:   PetscFunctionBegin;
408:   PetscAssertPointer(level, 2);
409:   *level = pc->kspnestlevel;
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: /*@
414:   PCCreate - Creates a preconditioner context, `PC`

416:   Collective

418:   Input Parameter:
419: . comm - MPI communicator

421:   Output Parameter:
422: . newpc - location to put the preconditioner context

424:   Level: developer

426:   Note:
427:   The default preconditioner for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC`
428:   in parallel. For dense matrices it is always `PCNONE`.

430: .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()`
431: @*/
432: PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
433: {
434:   PC pc;

436:   PetscFunctionBegin;
437:   PetscAssertPointer(newpc, 2);
438:   PetscCall(PCInitializePackage());

440:   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
441:   pc->mat                  = NULL;
442:   pc->pmat                 = NULL;
443:   pc->setupcalled          = 0;
444:   pc->setfromoptionscalled = 0;
445:   pc->data                 = NULL;
446:   pc->diagonalscale        = PETSC_FALSE;
447:   pc->diagonalscaleleft    = NULL;
448:   pc->diagonalscaleright   = NULL;

450:   pc->modifysubmatrices  = NULL;
451:   pc->modifysubmatricesP = NULL;

453:   *newpc = pc;
454:   PetscFunctionReturn(PETSC_SUCCESS);
455: }

457: /*@
458:   PCApply - Applies the preconditioner to a vector.

460:   Collective

462:   Input Parameters:
463: + pc - the preconditioner context
464: - x  - input vector

466:   Output Parameter:
467: . y - output vector

469:   Level: developer

471: .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
472: @*/
473: PetscErrorCode PCApply(PC pc, Vec x, Vec y)
474: {
475:   PetscInt m, n, mv, nv;

477:   PetscFunctionBegin;
481:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
482:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
483:   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
484:   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
485:   PetscCall(VecGetLocalSize(x, &mv));
486:   PetscCall(VecGetLocalSize(y, &nv));
487:   /* check pmat * y = x is feasible */
488:   PetscCheck(mv == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local rows %" PetscInt_FMT " does not equal input vector size %" PetscInt_FMT, m, mv);
489:   PetscCheck(nv == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local columns %" PetscInt_FMT " does not equal output vector size %" PetscInt_FMT, n, nv);
490:   PetscCall(VecSetErrorIfLocked(y, 3));

492:   PetscCall(PCSetUp(pc));
493:   PetscCall(VecLockReadPush(x));
494:   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
495:   PetscUseTypeMethod(pc, apply, x, y);
496:   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
497:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
498:   PetscCall(VecLockReadPop(x));
499:   PetscFunctionReturn(PETSC_SUCCESS);
500: }

502: /*@
503:   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.

505:   Collective

507:   Input Parameters:
508: + pc - the preconditioner context
509: - X  - block of input vectors

511:   Output Parameter:
512: . Y - block of output vectors

514:   Level: developer

516: .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
517: @*/
518: PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
519: {
520:   Mat       A;
521:   Vec       cy, cx;
522:   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
523:   PetscBool match;

525:   PetscFunctionBegin;
529:   PetscCheckSameComm(pc, 1, X, 2);
530:   PetscCheckSameComm(pc, 1, Y, 3);
531:   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
532:   PetscCall(PCGetOperators(pc, NULL, &A));
533:   PetscCall(MatGetLocalSize(A, &m3, &n3));
534:   PetscCall(MatGetLocalSize(X, &m2, &n2));
535:   PetscCall(MatGetLocalSize(Y, &m1, &n1));
536:   PetscCall(MatGetSize(A, &M3, &N3));
537:   PetscCall(MatGetSize(X, &M2, &N2));
538:   PetscCall(MatGetSize(Y, &M1, &N1));
539:   PetscCheck(n1 == n2 && N1 == N2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of input vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and block of output vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", n2, N2, n1, N1);
540:   PetscCheck(m2 == m3 && M2 == M3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of input vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m2, M2, m3, M3, n3, N3);
541:   PetscCheck(m1 == n3 && M1 == N3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of output vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m1, M1, m3, M3, n3, N3);
542:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
543:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
544:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
545:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
546:   PetscCall(PCSetUp(pc));
547:   if (pc->ops->matapply) {
548:     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
549:     PetscUseTypeMethod(pc, matapply, X, Y);
550:     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
551:   } else {
552:     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
553:     for (n1 = 0; n1 < N1; ++n1) {
554:       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
555:       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
556:       PetscCall(PCApply(pc, cx, cy));
557:       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
558:       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
559:     }
560:   }
561:   PetscFunctionReturn(PETSC_SUCCESS);
562: }

564: /*@
565:   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.

567:   Collective

569:   Input Parameters:
570: + pc - the preconditioner context
571: - x  - input vector

573:   Output Parameter:
574: . y - output vector

576:   Level: developer

578:   Note:
579:   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.

581: .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
582: @*/
583: PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
584: {
585:   PetscFunctionBegin;
589:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
590:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
591:   PetscCall(PCSetUp(pc));
592:   PetscCall(VecLockReadPush(x));
593:   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
594:   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
595:   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
596:   PetscCall(VecLockReadPop(x));
597:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
598:   PetscFunctionReturn(PETSC_SUCCESS);
599: }

601: /*@
602:   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.

604:   Collective

606:   Input Parameters:
607: + pc - the preconditioner context
608: - x  - input vector

610:   Output Parameter:
611: . y - output vector

613:   Level: developer

615:   Note:
616:   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.

618: .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
619: @*/
620: PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
621: {
622:   PetscFunctionBegin;
626:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
627:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
628:   PetscCall(PCSetUp(pc));
629:   PetscCall(VecLockReadPush(x));
630:   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
631:   PetscUseTypeMethod(pc, applysymmetricright, x, y);
632:   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
633:   PetscCall(VecLockReadPop(x));
634:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
635:   PetscFunctionReturn(PETSC_SUCCESS);
636: }

638: /*@
639:   PCApplyTranspose - Applies the transpose of preconditioner to a vector.

641:   Collective

643:   Input Parameters:
644: + pc - the preconditioner context
645: - x  - input vector

647:   Output Parameter:
648: . y - output vector

650:   Level: developer

652:   Note:
653:   For complex numbers this applies the non-Hermitian transpose.

655:   Developer Note:
656:   We need to implement a `PCApplyHermitianTranspose()`

658: .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
659: @*/
660: PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
661: {
662:   PetscFunctionBegin;
666:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
667:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
668:   PetscCall(PCSetUp(pc));
669:   PetscCall(VecLockReadPush(x));
670:   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
671:   PetscUseTypeMethod(pc, applytranspose, x, y);
672:   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
673:   PetscCall(VecLockReadPop(x));
674:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
675:   PetscFunctionReturn(PETSC_SUCCESS);
676: }

678: /*@
679:   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation

681:   Collective

683:   Input Parameter:
684: . pc - the preconditioner context

686:   Output Parameter:
687: . flg - `PETSC_TRUE` if a transpose operation is defined

689:   Level: developer

691: .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
692: @*/
693: PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
694: {
695:   PetscFunctionBegin;
697:   PetscAssertPointer(flg, 2);
698:   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
699:   else *flg = PETSC_FALSE;
700:   PetscFunctionReturn(PETSC_SUCCESS);
701: }

703: /*@
704:   PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.

706:   Collective

708:   Input Parameters:
709: + pc   - the preconditioner context
710: . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
711: . x    - input vector
712: - work - work vector

714:   Output Parameter:
715: . y - output vector

717:   Level: developer

719:   Note:
720:   If the `PC` has had `PCSetDiagonalScale()` set then $ D M A D^{-1} $ for left preconditioning or $ D A M D^{-1} $ is actually applied.
721:   The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.

723: .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
724: @*/
725: PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
726: {
727:   PetscFunctionBegin;
733:   PetscCheckSameComm(pc, 1, x, 3);
734:   PetscCheckSameComm(pc, 1, y, 4);
735:   PetscCheckSameComm(pc, 1, work, 5);
736:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
737:   PetscCheck(side == PC_LEFT || side == PC_SYMMETRIC || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right, left, or symmetric");
738:   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
739:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));

741:   PetscCall(PCSetUp(pc));
742:   if (pc->diagonalscale) {
743:     if (pc->ops->applyBA) {
744:       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
745:       PetscCall(VecDuplicate(x, &work2));
746:       PetscCall(PCDiagonalScaleRight(pc, x, work2));
747:       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
748:       PetscCall(PCDiagonalScaleLeft(pc, y, y));
749:       PetscCall(VecDestroy(&work2));
750:     } else if (side == PC_RIGHT) {
751:       PetscCall(PCDiagonalScaleRight(pc, x, y));
752:       PetscCall(PCApply(pc, y, work));
753:       PetscCall(MatMult(pc->mat, work, y));
754:       PetscCall(PCDiagonalScaleLeft(pc, y, y));
755:     } else if (side == PC_LEFT) {
756:       PetscCall(PCDiagonalScaleRight(pc, x, y));
757:       PetscCall(MatMult(pc->mat, y, work));
758:       PetscCall(PCApply(pc, work, y));
759:       PetscCall(PCDiagonalScaleLeft(pc, y, y));
760:     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
761:   } else {
762:     if (pc->ops->applyBA) {
763:       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
764:     } else if (side == PC_RIGHT) {
765:       PetscCall(PCApply(pc, x, work));
766:       PetscCall(MatMult(pc->mat, work, y));
767:     } else if (side == PC_LEFT) {
768:       PetscCall(MatMult(pc->mat, x, work));
769:       PetscCall(PCApply(pc, work, y));
770:     } else if (side == PC_SYMMETRIC) {
771:       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
772:       PetscCall(PCApplySymmetricRight(pc, x, work));
773:       PetscCall(MatMult(pc->mat, work, y));
774:       PetscCall(VecCopy(y, work));
775:       PetscCall(PCApplySymmetricLeft(pc, work, y));
776:     }
777:   }
778:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
779:   PetscFunctionReturn(PETSC_SUCCESS);
780: }

782: /*@
783:   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
784:   and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
785:   NOT $(B*A)^T = A^T*B^T$.

787:   Collective

789:   Input Parameters:
790: + pc   - the preconditioner context
791: . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
792: . x    - input vector
793: - work - work vector

795:   Output Parameter:
796: . y - output vector

798:   Level: developer

800:   Note:
801:   This routine is used internally so that the same Krylov code can be used to solve $A x = b$ and $A^T x = b$, with a preconditioner
802:   defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$

804: .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
805: @*/
806: PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
807: {
808:   PetscFunctionBegin;
813:   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
814:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
815:   if (pc->ops->applyBAtranspose) {
816:     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
817:     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
818:     PetscFunctionReturn(PETSC_SUCCESS);
819:   }
820:   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");

822:   PetscCall(PCSetUp(pc));
823:   if (side == PC_RIGHT) {
824:     PetscCall(PCApplyTranspose(pc, x, work));
825:     PetscCall(MatMultTranspose(pc->mat, work, y));
826:   } else if (side == PC_LEFT) {
827:     PetscCall(MatMultTranspose(pc->mat, x, work));
828:     PetscCall(PCApplyTranspose(pc, work, y));
829:   }
830:   /* add support for PC_SYMMETRIC */
831:   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
832:   PetscFunctionReturn(PETSC_SUCCESS);
833: }

835: /*@
836:   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
837:   built-in fast application of Richardson's method.

839:   Not Collective

841:   Input Parameter:
842: . pc - the preconditioner

844:   Output Parameter:
845: . exists - `PETSC_TRUE` or `PETSC_FALSE`

847:   Level: developer

849: .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()`
850: @*/
851: PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
852: {
853:   PetscFunctionBegin;
855:   PetscAssertPointer(exists, 2);
856:   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
857:   else *exists = PETSC_FALSE;
858:   PetscFunctionReturn(PETSC_SUCCESS);
859: }

861: /*@
862:   PCApplyRichardson - Applies several steps of Richardson iteration with
863:   the particular preconditioner. This routine is usually used by the
864:   Krylov solvers and not the application code directly.

866:   Collective

868:   Input Parameters:
869: + pc        - the preconditioner context
870: . b         - the right-hand side
871: . w         - one work vector
872: . rtol      - relative decrease in residual norm convergence criteria
873: . abstol    - absolute residual norm convergence criteria
874: . dtol      - divergence residual norm increase criteria
875: . its       - the number of iterations to apply.
876: - guesszero - if the input x contains nonzero initial guess

878:   Output Parameters:
879: + outits - number of iterations actually used (for SOR this always equals its)
880: . reason - the reason the apply terminated
881: - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`

883:   Level: developer

885:   Notes:
886:   Most preconditioners do not support this function. Use the command
887:   `PCApplyRichardsonExists()` to determine if one does.

889:   Except for the `PCMG` this routine ignores the convergence tolerances
890:   and always runs for the number of iterations

892: .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
893: @*/
894: PetscErrorCode PCApplyRichardson(PC pc, Vec b, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason)
895: {
896:   PetscFunctionBegin;
901:   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
902:   PetscCall(PCSetUp(pc));
903:   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }

907: /*@
908:   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail

910:   Logically Collective

912:   Input Parameters:
913: + pc     - the preconditioner context
914: - reason - the reason it failedx

916:   Level: advanced

918: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
919: @*/
920: PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
921: {
922:   PetscFunctionBegin;
924:   pc->failedreason = reason;
925:   PetscFunctionReturn(PETSC_SUCCESS);
926: }

928: /*@
929:   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail

931:   Logically Collective

933:   Input Parameter:
934: . pc - the preconditioner context

936:   Output Parameter:
937: . reason - the reason it failed

939:   Level: advanced

941:   Note:
942:   This is the maximum over reason over all ranks in the PC communicator. It is only valid after
943:   a call `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`.
944:   It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()`

946: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()`
947: @*/
948: PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
949: {
950:   PetscFunctionBegin;
952:   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
953:   else *reason = pc->failedreason;
954:   PetscFunctionReturn(PETSC_SUCCESS);
955: }

957: /*@
958:   PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank

960:   Not Collective

962:   Input Parameter:
963: . pc - the preconditioner context

965:   Output Parameter:
966: . reason - the reason it failed

968:   Level: advanced

970:   Note:
971:   Different processes may have different reasons or no reason, see `PCGetFailedReason()`

973: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()`
974: @*/
975: PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason)
976: {
977:   PetscFunctionBegin;
979:   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
980:   else *reason = pc->failedreason;
981:   PetscFunctionReturn(PETSC_SUCCESS);
982: }

984: /*@
985:   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`

987:   Collective

989:   Input Parameter:
990: . pc - the preconditioner context

992:   Level: advanced

994:   Note:
995:   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
996:   makes them have a common value (failure if any MPI process had a failure).

998: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`
999: @*/
1000: PetscErrorCode PCReduceFailedReason(PC pc)
1001: {
1002:   PetscInt buf;

1004:   PetscFunctionBegin;
1006:   buf = (PetscInt)pc->failedreason;
1007:   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
1008:   pc->failedreason = (PCFailedReason)buf;
1009:   PetscFunctionReturn(PETSC_SUCCESS);
1010: }

1012: /*  Next line needed to deactivate KSP_Solve logging */
1013: #include <petsc/private/kspimpl.h>

1015: /*
1016:       a setupcall of 0 indicates never setup,
1017:                      1 indicates has been previously setup
1018:                     -1 indicates a PCSetUp() was attempted and failed
1019: */
1020: /*@
1021:   PCSetUp - Prepares for the use of a preconditioner.

1023:   Collective

1025:   Input Parameter:
1026: . pc - the preconditioner context

1028:   Level: developer

1030: .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`
1031: @*/
1032: PetscErrorCode PCSetUp(PC pc)
1033: {
1034:   const char      *def;
1035:   PetscObjectState matstate, matnonzerostate;

1037:   PetscFunctionBegin;
1039:   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");

1041:   if (pc->setupcalled && pc->reusepreconditioner) {
1042:     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
1043:     PetscFunctionReturn(PETSC_SUCCESS);
1044:   }

1046:   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
1047:   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
1048:   if (!pc->setupcalled) {
1049:     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
1050:     pc->flag = DIFFERENT_NONZERO_PATTERN;
1051:   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
1052:   else {
1053:     if (matnonzerostate != pc->matnonzerostate) {
1054:       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
1055:       pc->flag = DIFFERENT_NONZERO_PATTERN;
1056:     } else {
1057:       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
1058:       pc->flag = SAME_NONZERO_PATTERN;
1059:     }
1060:   }
1061:   pc->matstate        = matstate;
1062:   pc->matnonzerostate = matnonzerostate;

1064:   if (!((PetscObject)pc)->type_name) {
1065:     PetscCall(PCGetDefaultType_Private(pc, &def));
1066:     PetscCall(PCSetType(pc, def));
1067:   }

1069:   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
1070:   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
1071:   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
1072:   if (pc->ops->setup) {
1073:     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
1074:     PetscCall(KSPInitializePackage());
1075:     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
1076:     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1077:     PetscUseTypeMethod(pc, setup);
1078:     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
1079:     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
1080:   }
1081:   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1082:   if (!pc->setupcalled) pc->setupcalled = 1;
1083:   PetscFunctionReturn(PETSC_SUCCESS);
1084: }

1086: /*@
1087:   PCSetUpOnBlocks - Sets up the preconditioner for each block in
1088:   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
1089:   methods.

1091:   Collective

1093:   Input Parameter:
1094: . pc - the preconditioner context

1096:   Level: developer

1098:   Note:
1099:   For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1100:   called on the outer `PC`, this routine ensures it is called.

1102: .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
1103: @*/
1104: PetscErrorCode PCSetUpOnBlocks(PC pc)
1105: {
1106:   PetscFunctionBegin;
1108:   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
1109:   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1110:   PetscUseTypeMethod(pc, setuponblocks);
1111:   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
1112:   PetscFunctionReturn(PETSC_SUCCESS);
1113: }

1115: /*@C
1116:   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
1117:   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`

1119:   Logically Collective

1121:   Input Parameters:
1122: + pc   - the preconditioner context
1123: . func - routine for modifying the submatrices
1124: - ctx  - optional user-defined context (may be `NULL`)

1126:   Calling sequence of `func`:
1127: + pc     - the preconditioner context
1128: . nsub   - number of index sets
1129: . row    - an array of index sets that contain the global row numbers
1130:          that comprise each local submatrix
1131: . col    - an array of index sets that contain the global column numbers
1132:          that comprise each local submatrix
1133: . submat - array of local submatrices
1134: - ctx    - optional user-defined context for private data for the
1135:          user-defined func routine (may be `NULL`)

1137:   Level: advanced

1139:   Notes:
1140:   The basic submatrices are extracted from the preconditioner matrix as
1141:   usual; the user can then alter these (for example, to set different boundary
1142:   conditions for each submatrix) before they are used for the local solves.

1144:   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1145:   `KSPSolve()`.

1147:   A routine set by `PCSetModifySubMatrices()` is currently called within
1148:   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
1149:   preconditioners.  All other preconditioners ignore this routine.

1151: .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
1152: @*/
1153: PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1154: {
1155:   PetscFunctionBegin;
1157:   pc->modifysubmatrices  = func;
1158:   pc->modifysubmatricesP = ctx;
1159:   PetscFunctionReturn(PETSC_SUCCESS);
1160: }

1162: /*@C
1163:   PCModifySubMatrices - Calls an optional user-defined routine within
1164:   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.

1166:   Collective

1168:   Input Parameters:
1169: + pc     - the preconditioner context
1170: . nsub   - the number of local submatrices
1171: . row    - an array of index sets that contain the global row numbers
1172:          that comprise each local submatrix
1173: . col    - an array of index sets that contain the global column numbers
1174:          that comprise each local submatrix
1175: . submat - array of local submatrices
1176: - ctx    - optional user-defined context for private data for the
1177:          user-defined routine (may be `NULL`)

1179:   Output Parameter:
1180: . submat - array of local submatrices (the entries of which may
1181:             have been modified)

1183:   Level: developer

1185:   Note:
1186:   The user should NOT generally call this routine, as it will
1187:   automatically be called within certain preconditioners.

1189: .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`
1190: @*/
1191: PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1192: {
1193:   PetscFunctionBegin;
1195:   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
1196:   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
1197:   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
1198:   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: /*@
1203:   PCSetOperators - Sets the matrix associated with the linear system and
1204:   a (possibly) different one associated with the preconditioner.

1206:   Logically Collective

1208:   Input Parameters:
1209: + pc   - the preconditioner context
1210: . Amat - the matrix that defines the linear system
1211: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

1213:   Level: intermediate

1215:   Notes:
1216:   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.

1218:   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1219:   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1220:   on it and then pass it back in in your call to `KSPSetOperators()`.

1222:   More Notes about Repeated Solution of Linear Systems:
1223:   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
1224:   to zero after a linear solve; the user is completely responsible for
1225:   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
1226:   zero all elements of a matrix.

1228: .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
1229:  @*/
1230: PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1231: {
1232:   PetscInt m1, n1, m2, n2;

1234:   PetscFunctionBegin;
1238:   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
1239:   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
1240:   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
1241:     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
1242:     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
1243:     PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Amat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1);
1244:     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
1245:     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
1246:     PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Pmat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1);
1247:   }

1249:   if (Pmat != pc->pmat) {
1250:     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1251:     pc->matnonzerostate = -1;
1252:     pc->matstate        = -1;
1253:   }

1255:   /* reference first in case the matrices are the same */
1256:   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
1257:   PetscCall(MatDestroy(&pc->mat));
1258:   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
1259:   PetscCall(MatDestroy(&pc->pmat));
1260:   pc->mat  = Amat;
1261:   pc->pmat = Pmat;
1262:   PetscFunctionReturn(PETSC_SUCCESS);
1263: }

1265: /*@
1266:   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.

1268:   Logically Collective

1270:   Input Parameters:
1271: + pc   - the preconditioner context
1272: - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner

1274:   Level: intermediate

1276:   Note:
1277:   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1278:   prevents this.

1280: .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
1281:  @*/
1282: PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1283: {
1284:   PetscFunctionBegin;
1287:   pc->reusepreconditioner = flag;
1288:   PetscFunctionReturn(PETSC_SUCCESS);
1289: }

1291: /*@
1292:   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.

1294:   Not Collective

1296:   Input Parameter:
1297: . pc - the preconditioner context

1299:   Output Parameter:
1300: . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner

1302:   Level: intermediate

1304: .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1305:  @*/
1306: PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1307: {
1308:   PetscFunctionBegin;
1310:   PetscAssertPointer(flag, 2);
1311:   *flag = pc->reusepreconditioner;
1312:   PetscFunctionReturn(PETSC_SUCCESS);
1313: }

1315: /*@
1316:   PCGetOperators - Gets the matrix associated with the linear system and
1317:   possibly a different one associated with the preconditioner.

1319:   Not Collective, though parallel `Mat`s are returned if `pc` is parallel

1321:   Input Parameter:
1322: . pc - the preconditioner context

1324:   Output Parameters:
1325: + Amat - the matrix defining the linear system
1326: - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.

1328:   Level: intermediate

1330:   Note:
1331:   Does not increase the reference count of the matrices, so you should not destroy them

1333:   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1334:   are created in `PC` and returned to the user. In this case, if both operators
1335:   mat and pmat are requested, two DIFFERENT operators will be returned. If
1336:   only one is requested both operators in the PC will be the same (i.e. as
1337:   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
1338:   The user must set the sizes of the returned matrices and their type etc just
1339:   as if the user created them with `MatCreate()`. For example,

1341: .vb
1342:          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1343:            set size, type, etc of Amat

1345:          MatCreate(comm,&mat);
1346:          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1347:          PetscObjectDereference((PetscObject)mat);
1348:            set size, type, etc of Amat
1349: .ve

1351:   and

1353: .vb
1354:          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1355:            set size, type, etc of Amat and Pmat

1357:          MatCreate(comm,&Amat);
1358:          MatCreate(comm,&Pmat);
1359:          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1360:          PetscObjectDereference((PetscObject)Amat);
1361:          PetscObjectDereference((PetscObject)Pmat);
1362:            set size, type, etc of Amat and Pmat
1363: .ve

1365:   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1366:   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1367:   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1368:   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1369:   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1370:   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1371:   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
1372:   it can be created for you?

1374: .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
1375: @*/
1376: PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1377: {
1378:   PetscFunctionBegin;
1380:   if (Amat) {
1381:     if (!pc->mat) {
1382:       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
1383:         pc->mat = pc->pmat;
1384:         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1385:       } else { /* both Amat and Pmat are empty */
1386:         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1387:         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1388:           pc->pmat = pc->mat;
1389:           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1390:         }
1391:       }
1392:     }
1393:     *Amat = pc->mat;
1394:   }
1395:   if (Pmat) {
1396:     if (!pc->pmat) {
1397:       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
1398:         pc->pmat = pc->mat;
1399:         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1400:       } else {
1401:         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1402:         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1403:           pc->mat = pc->pmat;
1404:           PetscCall(PetscObjectReference((PetscObject)pc->mat));
1405:         }
1406:       }
1407:     }
1408:     *Pmat = pc->pmat;
1409:   }
1410:   PetscFunctionReturn(PETSC_SUCCESS);
1411: }

1413: /*@
1414:   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1415:   possibly a different one associated with the preconditioner have been set in the `PC`.

1417:   Not Collective, though the results on all processes should be the same

1419:   Input Parameter:
1420: . pc - the preconditioner context

1422:   Output Parameters:
1423: + mat  - the matrix associated with the linear system was set
1424: - pmat - matrix associated with the preconditioner was set, usually the same

1426:   Level: intermediate

1428: .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1429: @*/
1430: PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1431: {
1432:   PetscFunctionBegin;
1434:   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1435:   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1436:   PetscFunctionReturn(PETSC_SUCCESS);
1437: }

1439: /*@
1440:   PCFactorGetMatrix - Gets the factored matrix from the
1441:   preconditioner context.  This routine is valid only for the `PCLU`,
1442:   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.

1444:   Not Collective though `mat` is parallel if `pc` is parallel

1446:   Input Parameter:
1447: . pc - the preconditioner context

1449:   Output Parameters:
1450: . mat - the factored matrix

1452:   Level: advanced

1454:   Note:
1455:   Does not increase the reference count for `mat` so DO NOT destroy it

1457: .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
1458: @*/
1459: PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1460: {
1461:   PetscFunctionBegin;
1463:   PetscAssertPointer(mat, 2);
1464:   PetscCall(PCFactorSetUpMatSolverType(pc));
1465:   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
1466:   PetscFunctionReturn(PETSC_SUCCESS);
1467: }

1469: /*@
1470:   PCSetOptionsPrefix - Sets the prefix used for searching for all
1471:   `PC` options in the database.

1473:   Logically Collective

1475:   Input Parameters:
1476: + pc     - the preconditioner context
1477: - prefix - the prefix string to prepend to all `PC` option requests

1479:   Note:
1480:   A hyphen (-) must NOT be given at the beginning of the prefix name.
1481:   The first character of all runtime options is AUTOMATICALLY the
1482:   hyphen.

1484:   Level: advanced

1486: .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
1487: @*/
1488: PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1489: {
1490:   PetscFunctionBegin;
1492:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
1493:   PetscFunctionReturn(PETSC_SUCCESS);
1494: }

1496: /*@
1497:   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1498:   `PC` options in the database.

1500:   Logically Collective

1502:   Input Parameters:
1503: + pc     - the preconditioner context
1504: - prefix - the prefix string to prepend to all `PC` option requests

1506:   Note:
1507:   A hyphen (-) must NOT be given at the beginning of the prefix name.
1508:   The first character of all runtime options is AUTOMATICALLY the
1509:   hyphen.

1511:   Level: advanced

1513: .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
1514: @*/
1515: PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1516: {
1517:   PetscFunctionBegin;
1519:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
1520:   PetscFunctionReturn(PETSC_SUCCESS);
1521: }

1523: /*@
1524:   PCGetOptionsPrefix - Gets the prefix used for searching for all
1525:   PC options in the database.

1527:   Not Collective

1529:   Input Parameter:
1530: . pc - the preconditioner context

1532:   Output Parameter:
1533: . prefix - pointer to the prefix string used, is returned

1535:   Level: advanced

1537:   Fortran Note:
1538:   The user should pass in a string `prefix` of
1539:   sufficient length to hold the prefix.

1541: .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
1542: @*/
1543: PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1544: {
1545:   PetscFunctionBegin;
1547:   PetscAssertPointer(prefix, 2);
1548:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
1549:   PetscFunctionReturn(PETSC_SUCCESS);
1550: }

1552: /*
1553:    Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
1554:   preconditioners including BDDC and Eisentat that transform the equations before applying
1555:   the Krylov methods
1556: */
1557: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1558: {
1559:   PetscFunctionBegin;
1561:   PetscAssertPointer(change, 2);
1562:   *change = PETSC_FALSE;
1563:   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
1564:   PetscFunctionReturn(PETSC_SUCCESS);
1565: }

1567: /*@
1568:   PCPreSolve - Optional pre-solve phase, intended for any
1569:   preconditioner-specific actions that must be performed before
1570:   the iterative solve itself.

1572:   Collective

1574:   Input Parameters:
1575: + pc  - the preconditioner context
1576: - ksp - the Krylov subspace context

1578:   Level: developer

1580:   Example Usage:
1581: .vb
1582:     PCPreSolve(pc,ksp);
1583:     KSPSolve(ksp,b,x);
1584:     PCPostSolve(pc,ksp);
1585: .ve

1587:   Notes:
1588:   The pre-solve phase is distinct from the `PCSetUp()` phase.

1590:   `KSPSolve()` calls this directly, so is rarely called by the user.

1592: .seealso: [](ch_ksp), `PC`, `PCPostSolve()`
1593: @*/
1594: PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1595: {
1596:   Vec x, rhs;

1598:   PetscFunctionBegin;
1601:   pc->presolvedone++;
1602:   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
1603:   PetscCall(KSPGetSolution(ksp, &x));
1604:   PetscCall(KSPGetRhs(ksp, &rhs));

1606:   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
1607:   else if (pc->presolve) PetscCall(pc->presolve(pc, ksp));
1608:   PetscFunctionReturn(PETSC_SUCCESS);
1609: }

1611: /*@C
1612:   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1613:   preconditioner-specific actions that must be performed before
1614:   the iterative solve itself.

1616:   Logically Collective

1618:   Input Parameters:
1619: + pc       - the preconditioner object
1620: - presolve - the function to call before the solve

1622:   Calling sequence of `presolve`:
1623: + pc  - the `PC` context
1624: - ksp - the `KSP` context

1626:   Level: developer

1628: .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()`
1629: @*/
1630: PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1631: {
1632:   PetscFunctionBegin;
1634:   pc->presolve = presolve;
1635:   PetscFunctionReturn(PETSC_SUCCESS);
1636: }

1638: /*@
1639:   PCPostSolve - Optional post-solve phase, intended for any
1640:   preconditioner-specific actions that must be performed after
1641:   the iterative solve itself.

1643:   Collective

1645:   Input Parameters:
1646: + pc  - the preconditioner context
1647: - ksp - the Krylov subspace context

1649:   Example Usage:
1650: .vb
1651:     PCPreSolve(pc,ksp);
1652:     KSPSolve(ksp,b,x);
1653:     PCPostSolve(pc,ksp);
1654: .ve

1656:   Level: developer

1658:   Note:
1659:   `KSPSolve()` calls this routine directly, so it is rarely called by the user.

1661: .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
1662: @*/
1663: PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1664: {
1665:   Vec x, rhs;

1667:   PetscFunctionBegin;
1670:   pc->presolvedone--;
1671:   PetscCall(KSPGetSolution(ksp, &x));
1672:   PetscCall(KSPGetRhs(ksp, &rhs));
1673:   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
1674:   PetscFunctionReturn(PETSC_SUCCESS);
1675: }

1677: /*@
1678:   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.

1680:   Collective

1682:   Input Parameters:
1683: + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1684:            some related function before a call to `PCLoad()`.
1685: - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`

1687:   Level: intermediate

1689:   Note:
1690:   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.

1692: .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
1693: @*/
1694: PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1695: {
1696:   PetscBool isbinary;
1697:   PetscInt  classid;
1698:   char      type[256];

1700:   PetscFunctionBegin;
1703:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1704:   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");

1706:   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
1707:   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
1708:   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
1709:   PetscCall(PCSetType(newdm, type));
1710:   PetscTryTypeMethod(newdm, load, viewer);
1711:   PetscFunctionReturn(PETSC_SUCCESS);
1712: }

1714: #include <petscdraw.h>
1715: #if defined(PETSC_HAVE_SAWS)
1716: #include <petscviewersaws.h>
1717: #endif

1719: /*@
1720:   PCViewFromOptions - View from the `PC` based on options in the options database

1722:   Collective

1724:   Input Parameters:
1725: + A    - the `PC` context
1726: . obj  - Optional object that provides the options prefix
1727: - name - command line option

1729:   Level: intermediate

1731: .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1732: @*/
1733: PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1734: {
1735:   PetscFunctionBegin;
1737:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1738:   PetscFunctionReturn(PETSC_SUCCESS);
1739: }

1741: /*@
1742:   PCView - Prints information about the `PC`

1744:   Collective

1746:   Input Parameters:
1747: + pc     - the `PC` context
1748: - viewer - optional visualization context

1750:   Level: developer

1752:   Notes:
1753:   The available visualization contexts include
1754: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1755: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
1756:   output where only the first processor opens
1757:   the file. All other processors send their
1758:   data to the first processor to print.

1760:   The user can open an alternative visualization contexts with
1761:   `PetscViewerASCIIOpen()` (output to a specified file).

1763: .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
1764: @*/
1765: PetscErrorCode PCView(PC pc, PetscViewer viewer)
1766: {
1767:   PCType            cstr;
1768:   PetscViewerFormat format;
1769:   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1770: #if defined(PETSC_HAVE_SAWS)
1771:   PetscBool issaws;
1772: #endif

1774:   PetscFunctionBegin;
1776:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
1778:   PetscCheckSameComm(pc, 1, viewer, 2);

1780:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1781:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1782:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1783:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1784: #if defined(PETSC_HAVE_SAWS)
1785:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
1786: #endif

1788:   if (iascii) {
1789:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
1790:     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
1791:     PetscCall(PetscViewerASCIIPushTab(viewer));
1792:     PetscTryTypeMethod(pc, view, viewer);
1793:     PetscCall(PetscViewerASCIIPopTab(viewer));
1794:     if (pc->mat) {
1795:       PetscCall(PetscViewerGetFormat(viewer, &format));
1796:       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1797:         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1798:         pop = PETSC_TRUE;
1799:       }
1800:       if (pc->pmat == pc->mat) {
1801:         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
1802:         PetscCall(PetscViewerASCIIPushTab(viewer));
1803:         PetscCall(MatView(pc->mat, viewer));
1804:         PetscCall(PetscViewerASCIIPopTab(viewer));
1805:       } else {
1806:         if (pc->pmat) {
1807:           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
1808:         } else {
1809:           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
1810:         }
1811:         PetscCall(PetscViewerASCIIPushTab(viewer));
1812:         PetscCall(MatView(pc->mat, viewer));
1813:         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
1814:         PetscCall(PetscViewerASCIIPopTab(viewer));
1815:       }
1816:       if (pop) PetscCall(PetscViewerPopFormat(viewer));
1817:     }
1818:   } else if (isstring) {
1819:     PetscCall(PCGetType(pc, &cstr));
1820:     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1821:     PetscTryTypeMethod(pc, view, viewer);
1822:     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1823:     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1824:   } else if (isbinary) {
1825:     PetscInt    classid = PC_FILE_CLASSID;
1826:     MPI_Comm    comm;
1827:     PetscMPIInt rank;
1828:     char        type[256];

1830:     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1831:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1832:     if (rank == 0) {
1833:       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
1834:       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
1835:       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
1836:     }
1837:     PetscTryTypeMethod(pc, view, viewer);
1838:   } else if (isdraw) {
1839:     PetscDraw draw;
1840:     char      str[25];
1841:     PetscReal x, y, bottom, h;
1842:     PetscInt  n;

1844:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1845:     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
1846:     if (pc->mat) {
1847:       PetscCall(MatGetSize(pc->mat, &n, NULL));
1848:       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
1849:     } else {
1850:       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
1851:     }
1852:     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
1853:     bottom = y - h;
1854:     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1855:     PetscTryTypeMethod(pc, view, viewer);
1856:     PetscCall(PetscDrawPopCurrentPoint(draw));
1857: #if defined(PETSC_HAVE_SAWS)
1858:   } else if (issaws) {
1859:     PetscMPIInt rank;

1861:     PetscCall(PetscObjectName((PetscObject)pc));
1862:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1863:     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
1864:     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1865:     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1866: #endif
1867:   }
1868:   PetscFunctionReturn(PETSC_SUCCESS);
1869: }

1871: /*@C
1872:   PCRegister -  Adds a method (`PCType`) to the preconditioner package.

1874:   Not collective. No Fortran Support

1876:   Input Parameters:
1877: + sname    - name of a new user-defined solver
1878: - function - routine to create method context

1880:   Example Usage:
1881: .vb
1882:    PCRegister("my_solver", MySolverCreate);
1883: .ve

1885:   Then, your solver can be chosen with the procedural interface via
1886: $     PCSetType(pc, "my_solver")
1887:   or at runtime via the option
1888: $     -pc_type my_solver

1890:   Level: advanced

1892:   Note:
1893:   `PCRegister()` may be called multiple times to add several user-defined preconditioners.

1895: .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`
1896: @*/
1897: PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1898: {
1899:   PetscFunctionBegin;
1900:   PetscCall(PCInitializePackage());
1901:   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
1902:   PetscFunctionReturn(PETSC_SUCCESS);
1903: }

1905: static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1906: {
1907:   PC pc;

1909:   PetscFunctionBegin;
1910:   PetscCall(MatShellGetContext(A, &pc));
1911:   PetscCall(PCApply(pc, X, Y));
1912:   PetscFunctionReturn(PETSC_SUCCESS);
1913: }

1915: /*@
1916:   PCComputeOperator - Computes the explicit preconditioned operator.

1918:   Collective

1920:   Input Parameters:
1921: + pc      - the preconditioner object
1922: - mattype - the `MatType` to be used for the operator

1924:   Output Parameter:
1925: . mat - the explicit preconditioned operator

1927:   Level: advanced

1929:   Note:
1930:   This computation is done by applying the operators to columns of the identity matrix.
1931:   This routine is costly in general, and is recommended for use only with relatively small systems.
1932:   Currently, this routine uses a dense matrix format when `mattype` == `NULL`

1934: .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
1935: @*/
1936: PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1937: {
1938:   PetscInt N, M, m, n;
1939:   Mat      A, Apc;

1941:   PetscFunctionBegin;
1943:   PetscAssertPointer(mat, 3);
1944:   PetscCall(PCGetOperators(pc, &A, NULL));
1945:   PetscCall(MatGetLocalSize(A, &m, &n));
1946:   PetscCall(MatGetSize(A, &M, &N));
1947:   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
1948:   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
1949:   PetscCall(MatComputeOperator(Apc, mattype, mat));
1950:   PetscCall(MatDestroy(&Apc));
1951:   PetscFunctionReturn(PETSC_SUCCESS);
1952: }

1954: /*@
1955:   PCSetCoordinates - sets the coordinates of all the nodes on the local process

1957:   Collective

1959:   Input Parameters:
1960: + pc     - the solver context
1961: . dim    - the dimension of the coordinates 1, 2, or 3
1962: . nloc   - the blocked size of the coordinates array
1963: - coords - the coordinates array

1965:   Level: intermediate

1967:   Note:
1968:   `coords` is an array of the dim coordinates for the nodes on
1969:   the local processor, of size `dim`*`nloc`.
1970:   If there are 108 equation on a processor
1971:   for a displacement finite element discretization of elasticity (so
1972:   that there are nloc = 36 = 108/3 nodes) then the array must have 108
1973:   double precision values (ie, 3 * 36).  These x y z coordinates
1974:   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1975:   ... , N-1.z ].

1977: .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
1978: @*/
1979: PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1980: {
1981:   PetscFunctionBegin;
1984:   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
1985:   PetscFunctionReturn(PETSC_SUCCESS);
1986: }

1988: /*@
1989:   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)

1991:   Logically Collective

1993:   Input Parameter:
1994: . pc - the precondition context

1996:   Output Parameters:
1997: + num_levels     - the number of levels
1998: - interpolations - the interpolation matrices (size of `num_levels`-1)

2000:   Level: advanced

2002:   Developer Note:
2003:   Why is this here instead of in `PCMG` etc?

2005: .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2006: @*/
2007: PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2008: {
2009:   PetscFunctionBegin;
2011:   PetscAssertPointer(num_levels, 2);
2012:   PetscAssertPointer(interpolations, 3);
2013:   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
2014:   PetscFunctionReturn(PETSC_SUCCESS);
2015: }

2017: /*@
2018:   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)

2020:   Logically Collective

2022:   Input Parameter:
2023: . pc - the precondition context

2025:   Output Parameters:
2026: + num_levels      - the number of levels
2027: - coarseOperators - the coarse operator matrices (size of `num_levels`-1)

2029:   Level: advanced

2031:   Developer Note:
2032:   Why is this here instead of in `PCMG` etc?

2034: .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2035: @*/
2036: PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2037: {
2038:   PetscFunctionBegin;
2040:   PetscAssertPointer(num_levels, 2);
2041:   PetscAssertPointer(coarseOperators, 3);
2042:   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
2043:   PetscFunctionReturn(PETSC_SUCCESS);
2044: }