1 2 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3 4 /*@ 5 SNESApplyNPC - Calls SNESSolve() on preconditioner for the SNES 6 7 Collective on snes 8 9 Input Parameters: 10 + snes - the SNES context 11 . x - input vector 12 - f - optional; the function evaluation on x 13 14 Output Parameter: 15 . y - function vector, as set by SNESSetFunction() 16 17 Note: 18 SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is 19 with SNESComuteJacobian(). 20 21 Level: developer 22 23 .seealso: `SNESGetNPC()`, `SNESSetNPC()`, `SNESComputeFunction()` 24 @*/ 25 PetscErrorCode SNESApplyNPC(SNES snes, Vec x, Vec f, Vec y) { 26 PetscFunctionBegin; 27 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 28 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 29 PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 30 PetscCheckSameComm(snes, 1, x, 2); 31 PetscCheckSameComm(snes, 1, y, 4); 32 PetscCall(VecValidValues(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 PetscFunctionReturn(0); 41 } 42 PetscFunctionReturn(0); 43 } 44 45 PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes, Vec X, Vec F) { 46 /* This is to be used as an argument to SNESMF -- NOT as a "function" */ 47 SNESConvergedReason reason; 48 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(0); 58 } 59 60 /*@ 61 SNESGetNPCFunction - Gets the current function value and its norm from a nonlinear preconditioner after `SNESSolve()` has been called on that `SNES` 62 63 Collective on snes 64 65 Input Parameter: 66 . snes - the `SNES` context 67 68 Output Parameters: 69 + F - function vector 70 - fnorm - the norm of F 71 72 Level: developer 73 74 .seealso: `SNESGetNPC()`, `SNESSetNPC()`, `SNESComputeFunction()`, `SNESApplyNPC()`, `SNESSolve()` 75 @*/ 76 PetscErrorCode SNESGetNPCFunction(SNES snes, Vec F, PetscReal *fnorm) { 77 PCSide npcside; 78 SNESFunctionType functype; 79 SNESNormSchedule normschedule; 80 Vec FPC, XPC; 81 82 PetscFunctionBegin; 83 if (snes->npc) { 84 PetscCall(SNESGetNPCSide(snes->npc, &npcside)); 85 PetscCall(SNESGetFunctionType(snes->npc, &functype)); 86 PetscCall(SNESGetNormSchedule(snes->npc, &normschedule)); 87 88 /* check if the function is valid based upon how the inner solver is preconditioned */ 89 if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) { 90 PetscCall(SNESGetFunction(snes->npc, &FPC, NULL, NULL)); 91 if (FPC) { 92 if (fnorm) PetscCall(VecNorm(FPC, NORM_2, fnorm)); 93 PetscCall(VecCopy(FPC, F)); 94 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlinear preconditioner has no function"); 95 } else { 96 PetscCall(SNESGetSolution(snes->npc, &XPC)); 97 if (XPC) { 98 PetscCall(SNESComputeFunction(snes->npc, XPC, F)); 99 if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm)); 100 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlinear preconditioner has no solution"); 101 } 102 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "No nonlinear preconditioner set"); 103 PetscFunctionReturn(0); 104 } 105