xref: /petsc/src/snes/interface/snespc.c (revision 1d3ed18c5acfa79332a5a7519adf0cfe659e9e83)
132f3f7c2SPeter Brune 
232f3f7c2SPeter Brune #include <petsc-private/snesimpl.h>      /*I "petscsnes.h"  I*/
332f3f7c2SPeter Brune #include <petscdmshell.h>
432f3f7c2SPeter Brune 
532f3f7c2SPeter Brune 
632f3f7c2SPeter Brune #undef __FUNCT__
732f3f7c2SPeter Brune #define __FUNCT__ "SNESApplyPC"
832f3f7c2SPeter Brune /*@
9ed07d7d7SPeter Brune    SNESApplyPC - Calls the function that has been set with SNESSetFunction().
1032f3f7c2SPeter Brune 
1132f3f7c2SPeter Brune    Collective on SNES
1232f3f7c2SPeter Brune 
1332f3f7c2SPeter Brune    Input Parameters:
1432f3f7c2SPeter Brune +  snes - the SNES context
1532f3f7c2SPeter Brune -  x - input vector
1632f3f7c2SPeter Brune 
1732f3f7c2SPeter Brune    Output Parameter:
1832f3f7c2SPeter Brune .  y - function vector, as set by SNESSetFunction()
1932f3f7c2SPeter Brune 
2032f3f7c2SPeter Brune    Notes:
2132f3f7c2SPeter Brune    SNESComputeFunction() should be called on X before SNESApplyPC() is called, as it is
2232f3f7c2SPeter Brune    with SNESComuteJacobian().
2332f3f7c2SPeter Brune 
2432f3f7c2SPeter Brune    Level: developer
2532f3f7c2SPeter Brune 
2632f3f7c2SPeter Brune .keywords: SNES, nonlinear, compute, function
2732f3f7c2SPeter Brune 
2832f3f7c2SPeter Brune .seealso: SNESGetPC(),SNESSetPC(),SNESComputeFunction()
2932f3f7c2SPeter Brune @*/
3032f3f7c2SPeter Brune PetscErrorCode  SNESApplyPC(SNES snes,Vec x,Vec f,PetscReal *fnorm,Vec y)
3132f3f7c2SPeter Brune {
3232f3f7c2SPeter Brune   PetscErrorCode ierr;
3332f3f7c2SPeter Brune 
3432f3f7c2SPeter Brune   PetscFunctionBegin;
3532f3f7c2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3632f3f7c2SPeter Brune   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3732f3f7c2SPeter Brune   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
3832f3f7c2SPeter Brune   PetscCheckSameComm(snes,1,x,2);
3932f3f7c2SPeter Brune   PetscCheckSameComm(snes,1,y,3);
4032f3f7c2SPeter Brune   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
4132f3f7c2SPeter Brune   if (snes->pc) {
4232f3f7c2SPeter Brune     if (f) {
4332f3f7c2SPeter Brune       ierr = SNESSetInitialFunction(snes->pc,f);CHKERRQ(ierr);
4432f3f7c2SPeter Brune     }
4532f3f7c2SPeter Brune     if (fnorm) {
4632f3f7c2SPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc,*fnorm);CHKERRQ(ierr);
4732f3f7c2SPeter Brune     }
4832f3f7c2SPeter Brune     ierr = VecCopy(x,y);CHKERRQ(ierr);
4932f3f7c2SPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
5032f3f7c2SPeter Brune     ierr = SNESSolve(snes->pc,snes->vec_rhs,y);CHKERRQ(ierr);
5132f3f7c2SPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr);
5232f3f7c2SPeter Brune     ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr);
5332f3f7c2SPeter Brune     ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
5432f3f7c2SPeter Brune     PetscFunctionReturn(0);
5532f3f7c2SPeter Brune   }
5632f3f7c2SPeter Brune   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
5732f3f7c2SPeter Brune   PetscFunctionReturn(0);
5832f3f7c2SPeter Brune }
5932f3f7c2SPeter Brune 
6032f3f7c2SPeter Brune #undef __FUNCT__
6132f3f7c2SPeter Brune #define __FUNCT__ "SNESComputeFunctionDefaultPC"
62ed07d7d7SPeter Brune PetscErrorCode SNESComputeFunctionDefaultPC(SNES snes,Vec X,Vec F) {
6332f3f7c2SPeter Brune /* This is to be used as an argument to SNESMF -- NOT as a "function" */
64b7281c8aSPeter Brune   SNESConvergedReason reason;
6532f3f7c2SPeter Brune   PetscErrorCode      ierr;
66b7281c8aSPeter Brune 
6732f3f7c2SPeter Brune   PetscFunctionBegin;
68b7281c8aSPeter Brune   if (snes->pc) {
696950c336SPeter Brune     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
70b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
71b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
72b7281c8aSPeter Brune       ierr = SNESSetFunctionDomainError(snes);CHKERRQ(ierr);
73b7281c8aSPeter Brune     }
74b7281c8aSPeter Brune   } else {
75b7281c8aSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
76b7281c8aSPeter Brune   }
7732f3f7c2SPeter Brune PetscFunctionReturn(0);
7832f3f7c2SPeter Brune }
79ddd40ce5SPeter Brune 
80ddd40ce5SPeter Brune #undef __FUNCT__
81ddd40ce5SPeter Brune #define __FUNCT__ "SNESGetPCFunction"
82ddd40ce5SPeter Brune /*@
83ddd40ce5SPeter Brune    SNESGetPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.
84ddd40ce5SPeter Brune 
85ddd40ce5SPeter Brune    Collective on SNES
86ddd40ce5SPeter Brune 
87ddd40ce5SPeter Brune    Input Parameters:
88ddd40ce5SPeter Brune .  snes - the SNES context
89ddd40ce5SPeter Brune 
90ddd40ce5SPeter Brune    Output Parameter:
91ddd40ce5SPeter Brune .  F - function vector
92ddd40ce5SPeter Brune .  fnorm - the norm of F
93ddd40ce5SPeter Brune 
94ddd40ce5SPeter Brune    Level: developer
95ddd40ce5SPeter Brune 
96ddd40ce5SPeter Brune .keywords: SNES, nonlinear, function
97ddd40ce5SPeter Brune 
98ddd40ce5SPeter Brune .seealso: SNESGetPC(),SNESSetPC(),SNESComputeFunction(),SNESApplyPC(),SNESSolve()
99ddd40ce5SPeter Brune @*/
100ddd40ce5SPeter Brune PetscErrorCode SNESGetPCFunction(SNES snes,Vec F,PetscReal *fnorm)
101ddd40ce5SPeter Brune {
102ddd40ce5SPeter Brune   PetscErrorCode   ierr;
103ddd40ce5SPeter Brune   PCSide           npcside;
104ddd40ce5SPeter Brune   SNESFunctionType functype;
105ddd40ce5SPeter Brune   SNESNormSchedule normschedule;
106ddd40ce5SPeter Brune   Vec              FPC,XPC;
107ddd40ce5SPeter Brune 
108ddd40ce5SPeter Brune   PetscFunctionBegin;
109ddd40ce5SPeter Brune   if (snes->pc) {
110ddd40ce5SPeter Brune     ierr = SNESGetPCSide(snes->pc,&npcside);CHKERRQ(ierr);
111ddd40ce5SPeter Brune     ierr = SNESGetFunctionType(snes->pc,&functype);CHKERRQ(ierr);
112ddd40ce5SPeter Brune     ierr = SNESGetNormSchedule(snes->pc,&normschedule);CHKERRQ(ierr);
113ddd40ce5SPeter Brune 
114ddd40ce5SPeter Brune     /* check if the function is valid based upon how the inner solver is preconditioned */
115ddd40ce5SPeter Brune     if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) {
116ddd40ce5SPeter Brune       ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
117*1d3ed18cSPeter Brune       if (FPC) {
118ddd40ce5SPeter Brune         if (fnorm) {ierr = SNESGetFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);}
119ddd40ce5SPeter Brune         ierr = VecCopy(FPC,F);CHKERRQ(ierr);
120ddd40ce5SPeter Brune       } else {
121ddd40ce5SPeter Brune         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
122ddd40ce5SPeter Brune       }
123ddd40ce5SPeter Brune     } else {
124ddd40ce5SPeter Brune       ierr = SNESGetSolution(snes->pc,&XPC);CHKERRQ(ierr);
125ddd40ce5SPeter Brune       if (XPC) {
126ddd40ce5SPeter Brune         ierr = SNESComputeFunction(snes->pc,XPC,F);CHKERRQ(ierr);
127ddd40ce5SPeter Brune         if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
128ddd40ce5SPeter Brune       } else {
129ddd40ce5SPeter Brune         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution");
130ddd40ce5SPeter Brune       }
131ddd40ce5SPeter Brune     }
132ddd40ce5SPeter Brune   } else {
133ddd40ce5SPeter Brune     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set");
134ddd40ce5SPeter Brune   }
135ddd40ce5SPeter Brune 
136ddd40ce5SPeter Brune   PetscFunctionReturn(0);
137ddd40ce5SPeter Brune }
138