Actual source code: snespc.c
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: .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
24: @*/
25: PetscErrorCode SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
26: {
35: VecValidValues(x,2,PETSC_TRUE);
36: if (snes->npc) {
37: if (f) {
38: SNESSetInitialFunction(snes->npc,f);
39: }
40: VecCopy(x,y);
41: PetscLogEventBegin(SNES_NPCSolve,snes->npc,x,y,0);
42: SNESSolve(snes->npc,snes->vec_rhs,y);
43: PetscLogEventEnd(SNES_NPCSolve,snes->npc,x,y,0);
44: VecAYPX(y,-1.0,x);
45: return(0);
46: }
47: return(0);
48: }
50: PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes,Vec X,Vec F)
51: {
52: /* This is to be used as an argument to SNESMF -- NOT as a "function" */
53: SNESConvergedReason reason;
54: PetscErrorCode ierr;
57: if (snes->npc) {
58: SNESApplyNPC(snes,X,NULL,F);
59: SNESGetConvergedReason(snes->npc,&reason);
60: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
61: SNESSetFunctionDomainError(snes);
62: }
63: } else {
64: SNESComputeFunction(snes,X,F);
65: }
66: return(0);
67: }
69: /*@
70: SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.
72: Collective on SNES
74: Input Parameters:
75: . snes - the SNES context
77: Output Parameter:
78: + F - function vector
79: - fnorm - the norm of F
81: Level: developer
83: .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction(),SNESApplyNPC(),SNESSolve()
84: @*/
85: PetscErrorCode SNESGetNPCFunction(SNES snes,Vec F,PetscReal *fnorm)
86: {
87: PetscErrorCode ierr;
88: PCSide npcside;
89: SNESFunctionType functype;
90: SNESNormSchedule normschedule;
91: Vec FPC,XPC;
94: if (snes->npc) {
95: SNESGetNPCSide(snes->npc,&npcside);
96: SNESGetFunctionType(snes->npc,&functype);
97: SNESGetNormSchedule(snes->npc,&normschedule);
99: /* check if the function is valid based upon how the inner solver is preconditioned */
100: if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) {
101: SNESGetFunction(snes->npc,&FPC,NULL,NULL);
102: if (FPC) {
103: if (fnorm) {VecNorm(FPC,NORM_2,fnorm);}
104: VecCopy(FPC,F);
105: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
106: } else {
107: SNESGetSolution(snes->npc,&XPC);
108: if (XPC) {
109: SNESComputeFunction(snes->npc,XPC,F);
110: if (fnorm) {VecNorm(F,NORM_2,fnorm);}
111: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution");
112: }
113: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set");
114: return(0);
115: }