Actual source code: snespc.c

  1: #include <petsc/private/snesimpl.h>

  3: /*@
  4:   SNESApplyNPC - Calls `SNESSolve()` on the preconditioner for the `SNES`

  6:   Collective

  8:   Input Parameters:
  9: + snes - the `SNES` context
 10: . x    - input vector
 11: - f    - optional; the function evaluation on `x`

 13:   Output Parameter:
 14: . y - function vector, as set by `SNESSetFunction()`

 16:   Level: developer

 18:   Note:
 19:   `SNESComputeFunction()` should be called on `x` before `SNESApplyNPC()` is called, as it is
 20:   with `SNESComuteJacobian()`.

 22: .seealso: [](ch_snes), `SNES`, `SNESGetNPC()`, `SNESSetNPC()`, `SNESComputeFunction()`
 23: @*/
 24: PetscErrorCode SNESApplyNPC(SNES snes, Vec x, Vec f, Vec y)
 25: {
 26:   PetscFunctionBegin;
 30:   PetscCheckSameComm(snes, 1, x, 2);
 31:   PetscCheckSameComm(snes, 1, y, 4);
 32:   PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
 33:   if (snes->npc) {
 34:     if (f) PetscCall(SNESSetInitialFunction(snes->npc, f));
 35:     PetscCall(VecCopy(x, y));
 36:     PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, x, y, 0));
 37:     PetscCall(SNESSolve(snes->npc, snes->vec_rhs, y));
 38:     PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, x, y, 0));
 39:     PetscCall(VecAYPX(y, -1.0, x));
 40:   }
 41:   PetscFunctionReturn(PETSC_SUCCESS);
 42: }

 44: PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes, Vec X, Vec F)
 45: {
 46:   /* This is to be used as an argument to SNESMF -- NOT as a "function" */
 47:   SNESConvergedReason reason;

 49:   PetscFunctionBegin;
 50:   if (snes->npc) {
 51:     PetscCall(SNESApplyNPC(snes, X, NULL, F));
 52:     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
 53:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) PetscCall(SNESSetFunctionDomainError(snes));
 54:   } else {
 55:     PetscCall(SNESComputeFunction(snes, X, F));
 56:   }
 57:   PetscFunctionReturn(PETSC_SUCCESS);
 58: }

 60: /*@
 61:   SNESGetNPCFunction - Gets the current function value and its norm from a nonlinear preconditioner after `SNESSolve()` has been called on that `SNES`

 63:   Collective

 65:   Input Parameter:
 66: . snes - the `SNES` context

 68:   Output Parameters:
 69: + F     - function vector
 70: - fnorm - the norm of `F`

 72:   Level: developer

 74: .seealso: [](ch_snes), `SNES`, `SNESGetNPC()`, `SNESSetNPC()`, `SNESComputeFunction()`, `SNESApplyNPC()`, `SNESSolve()`
 75: @*/
 76: PetscErrorCode SNESGetNPCFunction(SNES snes, Vec F, PetscReal *fnorm)
 77: {
 78:   PCSide           npcside;
 79:   SNESFunctionType functype;
 80:   SNESNormSchedule normschedule;
 81:   Vec              FPC, XPC;

 83:   PetscFunctionBegin;
 85:   if (fnorm) PetscAssertPointer(fnorm, 3);
 86:   PetscCheck(snes->npc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "No nonlinear preconditioner set");
 87:   PetscCall(SNESGetNPCSide(snes->npc, &npcside));
 88:   PetscCall(SNESGetFunctionType(snes->npc, &functype));
 89:   PetscCall(SNESGetNormSchedule(snes->npc, &normschedule));

 91:   /* check if the function is valid based upon how the inner solver is preconditioned */
 92:   if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) {
 93:     PetscCall(SNESGetFunction(snes->npc, &FPC, NULL, NULL));
 94:     PetscCheck(FPC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlinear preconditioner has no function");
 95:     if (fnorm) PetscCall(VecNorm(FPC, NORM_2, fnorm));
 96:     PetscCall(VecCopy(FPC, F));
 97:   } else {
 98:     PetscCall(SNESGetSolution(snes->npc, &XPC));
 99:     PetscCheck(XPC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlinear preconditioner has no solution");
100:     PetscCall(SNESComputeFunction(snes->npc, XPC, F));
101:     if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm));
102:   }
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }