xref: /petsc/src/snes/interface/snespc.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
1 
2 #include <petsc/private/snesimpl.h>      /*I "petscsnes.h"  I*/
3 #include <petscdmshell.h>
4 
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "SNESApplyNPC"
8 /*@
9    SNESApplyNPC - Calls SNESSolve() on preconditioner for the SNES
10 
11    Collective on SNES
12 
13    Input Parameters:
14 +  snes - the SNES context
15 .  x - input vector
16 -  f - optional; the function evaluation on x
17 
18    Output Parameter:
19 .  y - function vector, as set by SNESSetFunction()
20 
21    Notes:
22    SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is
23    with SNESComuteJacobian().
24 
25    Level: developer
26 
27 .keywords: SNES, nonlinear, compute, function
28 
29 .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
30 @*/
31 PetscErrorCode  SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
32 {
33   PetscErrorCode ierr;
34 
35   PetscFunctionBegin;
36   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
37   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
38   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
39   PetscCheckSameComm(snes,1,x,2);
40   PetscCheckSameComm(snes,1,y,3);
41   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
42   if (snes->pc) {
43     if (f) {
44       ierr = SNESSetInitialFunction(snes->pc,f);CHKERRQ(ierr);
45     }
46     ierr = VecCopy(x,y);CHKERRQ(ierr);
47     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
48     ierr = SNESSolve(snes->pc,snes->vec_rhs,y);CHKERRQ(ierr);
49     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
50     ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr);
51     PetscFunctionReturn(0);
52   }
53   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "SNESComputeFunctionDefaultNPC"
58 PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes,Vec X,Vec F)
59 {
60 /* This is to be used as an argument to SNESMF -- NOT as a "function" */
61   SNESConvergedReason reason;
62   PetscErrorCode      ierr;
63 
64   PetscFunctionBegin;
65   if (snes->pc) {
66     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
67     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
68     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
69       ierr = SNESSetFunctionDomainError(snes);CHKERRQ(ierr);
70     }
71   } else {
72     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
73   }
74   PetscFunctionReturn(0);
75 }
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "SNESGetNPCFunction"
79 /*@
80    SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.
81 
82    Collective on SNES
83 
84    Input Parameters:
85 .  snes - the SNES context
86 
87    Output Parameter:
88 .  F - function vector
89 .  fnorm - the norm of F
90 
91    Level: developer
92 
93 .keywords: SNES, nonlinear, function
94 
95 .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction(),SNESApplyNPC(),SNESSolve()
96 @*/
97 PetscErrorCode SNESGetNPCFunction(SNES snes,Vec F,PetscReal *fnorm)
98 {
99   PetscErrorCode   ierr;
100   PCSide           npcside;
101   SNESFunctionType functype;
102   SNESNormSchedule normschedule;
103   Vec              FPC,XPC;
104 
105   PetscFunctionBegin;
106   if (snes->pc) {
107     ierr = SNESGetNPCSide(snes->pc,&npcside);CHKERRQ(ierr);
108     ierr = SNESGetFunctionType(snes->pc,&functype);CHKERRQ(ierr);
109     ierr = SNESGetNormSchedule(snes->pc,&normschedule);CHKERRQ(ierr);
110 
111     /* check if the function is valid based upon how the inner solver is preconditioned */
112     if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) {
113       ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
114       if (FPC) {
115         if (fnorm) {ierr = VecNorm(FPC,NORM_2,fnorm);CHKERRQ(ierr);}
116         ierr = VecCopy(FPC,F);CHKERRQ(ierr);
117       } else {
118         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
119       }
120     } else {
121       ierr = SNESGetSolution(snes->pc,&XPC);CHKERRQ(ierr);
122       if (XPC) {
123         ierr = SNESComputeFunction(snes->pc,XPC,F);CHKERRQ(ierr);
124         if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
125       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution");
126     }
127   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set");
128   PetscFunctionReturn(0);
129 }
130