xref: /petsc/src/snes/interface/snespc.c (revision bfc4c25e64040512a7cdc7bcb995f785926f51e6)
1 
2 #include <petsc/private/snesimpl.h>      /*I "petscsnes.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESApplyNPC"
6 /*@
7    SNESApplyNPC - Calls SNESSolve() on preconditioner for the SNES
8 
9    Collective on SNES
10 
11    Input Parameters:
12 +  snes - the SNES context
13 .  x - input vector
14 -  f - optional; the function evaluation on x
15 
16    Output Parameter:
17 .  y - function vector, as set by SNESSetFunction()
18 
19    Notes:
20    SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is
21    with SNESComuteJacobian().
22 
23    Level: developer
24 
25 .keywords: SNES, nonlinear, compute, function
26 
27 .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
28 @*/
29 PetscErrorCode  SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
30 {
31   PetscErrorCode ierr;
32 
33   PetscFunctionBegin;
34   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
35   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
36   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
37   PetscCheckSameComm(snes,1,x,2);
38   PetscCheckSameComm(snes,1,y,3);
39   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
40   if (snes->pc) {
41     if (f) {
42       ierr = SNESSetInitialFunction(snes->pc,f);CHKERRQ(ierr);
43     }
44     ierr = VecCopy(x,y);CHKERRQ(ierr);
45     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
46     ierr = SNESSolve(snes->pc,snes->vec_rhs,y);CHKERRQ(ierr);
47     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
48     ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr);
49     PetscFunctionReturn(0);
50   }
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "SNESComputeFunctionDefaultNPC"
56 PetscErrorCode SNESComputeFunctionDefaultNPC(SNES snes,Vec X,Vec F)
57 {
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 = SNESApplyNPC(snes,X,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__ "SNESGetNPCFunction"
77 /*@
78    SNESGetNPCFunction - 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: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction(),SNESApplyNPC(),SNESSolve()
94 @*/
95 PetscErrorCode SNESGetNPCFunction(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 = SNESGetNPCSide(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 = VecNorm(FPC,NORM_2,fnorm);CHKERRQ(ierr);}
114         ierr = VecCopy(FPC,F);CHKERRQ(ierr);
115       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
116     } else {
117       ierr = SNESGetSolution(snes->pc,&XPC);CHKERRQ(ierr);
118       if (XPC) {
119         ierr = SNESComputeFunction(snes->pc,XPC,F);CHKERRQ(ierr);
120         if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
121       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution");
122     }
123   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set");
124   PetscFunctionReturn(0);
125 }
126