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