xref: /petsc/src/snes/interface/snespc.c (revision feff33ee0b5b037fa8f9f294dede656a2f85cc47)
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    Notes:
18    SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is
19    with SNESComuteJacobian().
20 
21    Level: developer
22 
23 .keywords: SNES, nonlinear, compute, function
24 
25 .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
26 @*/
27 PetscErrorCode  SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
28 {
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
33   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
34   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
35   PetscCheckSameComm(snes,1,x,2);
36   PetscCheckSameComm(snes,1,y,3);
37   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
38   if (snes->npc) {
39     if (f) {
40       ierr = SNESSetInitialFunction(snes->npc,f);CHKERRQ(ierr);
41     }
42     ierr = VecCopy(x,y);CHKERRQ(ierr);
43     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,x,y,0);CHKERRQ(ierr);
44     ierr = SNESSolve(snes->npc,snes->vec_rhs,y);CHKERRQ(ierr);
45     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,x,y,0);CHKERRQ(ierr);
46     ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr);
47     PetscFunctionReturn(0);
48   }
49   PetscFunctionReturn(0);
50 }
51 
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;
57 
58   PetscFunctionBegin;
59   if (snes->npc) {
60     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
61     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
62     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
63       ierr = SNESSetFunctionDomainError(snes);CHKERRQ(ierr);
64     }
65   } else {
66     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
67   }
68   PetscFunctionReturn(0);
69 }
70 
71 /*@
72    SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.
73 
74    Collective on SNES
75 
76    Input Parameters:
77 .  snes - the SNES context
78 
79    Output Parameter:
80 .  F - function vector
81 .  fnorm - the norm of F
82 
83    Level: developer
84 
85 .keywords: SNES, nonlinear, function
86 
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;
96 
97   PetscFunctionBegin;
98   if (snes->npc) {
99     ierr = SNESGetNPCSide(snes->npc,&npcside);CHKERRQ(ierr);
100     ierr = SNESGetFunctionType(snes->npc,&functype);CHKERRQ(ierr);
101     ierr = SNESGetNormSchedule(snes->npc,&normschedule);CHKERRQ(ierr);
102 
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       ierr = SNESGetFunction(snes->npc,&FPC,NULL,NULL);CHKERRQ(ierr);
106       if (FPC) {
107         if (fnorm) {ierr = VecNorm(FPC,NORM_2,fnorm);CHKERRQ(ierr);}
108         ierr = VecCopy(FPC,F);CHKERRQ(ierr);
109       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
110     } else {
111       ierr = SNESGetSolution(snes->npc,&XPC);CHKERRQ(ierr);
112       if (XPC) {
113         ierr = SNESComputeFunction(snes->npc,XPC,F);CHKERRQ(ierr);
114         if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
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   PetscFunctionReturn(0);
119 }
120