Actual source code: snespc.c
petsc-3.11.4 2019-09-28
2: #include <petsc/private/snesimpl.h>
4: /*@
5: SNESApplyNPC - Calls SNESSolve() on preconditioner for the SNES
7: Collective on SNES
9: Input Parameters:
10: + snes - the SNES context
11: . x - input vector
12: - f - optional; the function evaluation on x
14: Output Parameter:
15: . y - function vector, as set by SNESSetFunction()
17: Notes:
18: SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is
19: with SNESComuteJacobian().
21: Level: developer
23: .keywords: SNES, nonlinear, compute, function
25: .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
26: @*/
27: PetscErrorCode SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
28: {
37: VecValidValues(x,2,PETSC_TRUE);
38: if (snes->npc) {
39: if (f) {
40: SNESSetInitialFunction(snes->npc,f);
41: }
42: VecCopy(x,y);
43: PetscLogEventBegin(SNES_NPCSolve,snes->npc,x,y,0);
44: SNESSolve(snes->npc,snes->vec_rhs,y);
45: PetscLogEventEnd(SNES_NPCSolve,snes->npc,x,y,0);
46: VecAYPX(y,-1.0,x);
47: return(0);
48: }
49: return(0);
50: }
52: PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes,Vec X,Vec F)
53: {
54: /* This is to be used as an argument to SNESMF -- NOT as a "function" */
55: SNESConvergedReason reason;
56: PetscErrorCode ierr;
59: if (snes->npc) {
60: SNESApplyNPC(snes,X,NULL,F);
61: SNESGetConvergedReason(snes->npc,&reason);
62: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
63: SNESSetFunctionDomainError(snes);
64: }
65: } else {
66: SNESComputeFunction(snes,X,F);
67: }
68: return(0);
69: }
71: /*@
72: SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.
74: Collective on SNES
76: Input Parameters:
77: . snes - the SNES context
79: Output Parameter:
80: . F - function vector
81: . fnorm - the norm of F
83: Level: developer
85: .keywords: SNES, nonlinear, function
87: .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction(),SNESApplyNPC(),SNESSolve()
88: @*/
89: PetscErrorCode SNESGetNPCFunction(SNES snes,Vec F,PetscReal *fnorm)
90: {
91: PetscErrorCode ierr;
92: PCSide npcside;
93: SNESFunctionType functype;
94: SNESNormSchedule normschedule;
95: Vec FPC,XPC;
98: if (snes->npc) {
99: SNESGetNPCSide(snes->npc,&npcside);
100: SNESGetFunctionType(snes->npc,&functype);
101: SNESGetNormSchedule(snes->npc,&normschedule);
103: /* check if the function is valid based upon how the inner solver is preconditioned */
104: if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) {
105: SNESGetFunction(snes->npc,&FPC,NULL,NULL);
106: if (FPC) {
107: if (fnorm) {VecNorm(FPC,NORM_2,fnorm);}
108: VecCopy(FPC,F);
109: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
110: } else {
111: SNESGetSolution(snes->npc,&XPC);
112: if (XPC) {
113: SNESComputeFunction(snes->npc,XPC,F);
114: if (fnorm) {VecNorm(F,NORM_2,fnorm);}
115: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution");
116: }
117: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set");
118: return(0);
119: }