Actual source code: snespc.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/snesimpl.h>      /*I "petscsnes.h"  I*/
  3: #include <petscdmshell.h>


  8: /*@
  9:    SNESApplyNPC - Calls SNESSolve() on preconditioner for the SNES

 11:    Collective on SNES

 13:    Input Parameters:
 14: +  snes - the SNES context
 15: .  x - input vector
 16: -  f - optional; the function evaluation on x

 18:    Output Parameter:
 19: .  y - function vector, as set by SNESSetFunction()

 21:    Notes:
 22:    SNESComputeFunction() should be called on x before SNESApplyNPC() is called, as it is
 23:    with SNESComuteJacobian().

 25:    Level: developer

 27: .keywords: SNES, nonlinear, compute, function

 29: .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction()
 30: @*/
 31: PetscErrorCode  SNESApplyNPC(SNES snes,Vec x,Vec f,Vec y)
 32: {

 41:   VecValidValues(x,2,PETSC_TRUE);
 42:   if (snes->pc) {
 43:     if (f) {
 44:       SNESSetInitialFunction(snes->pc,f);
 45:     }
 46:     VecCopy(x,y);
 47:     PetscLogEventBegin(SNES_NPCSolve,snes->pc,x,y,0);
 48:     SNESSolve(snes->pc,snes->vec_rhs,y);
 49:     PetscLogEventEnd(SNES_NPCSolve,snes->pc,x,y,0);
 50:     VecAYPX(y,-1.0,x);
 51:     return(0);
 52:   }
 53:   return(0);
 54: }

 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;

 65:   if (snes->pc) {
 66:     SNESApplyNPC(snes,X,NULL,F);
 67:     SNESGetConvergedReason(snes->pc,&reason);
 68:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
 69:       SNESSetFunctionDomainError(snes);
 70:     }
 71:   } else {
 72:     SNESComputeFunction(snes,X,F);
 73:   }
 74:   return(0);
 75: }

 79: /*@
 80:    SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called.

 82:    Collective on SNES

 84:    Input Parameters:
 85: .  snes - the SNES context

 87:    Output Parameter:
 88: .  F - function vector
 89: .  fnorm - the norm of F

 91:    Level: developer

 93: .keywords: SNES, nonlinear, function

 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;

106:   if (snes->pc) {
107:     SNESGetNPCSide(snes->pc,&npcside);
108:     SNESGetFunctionType(snes->pc,&functype);
109:     SNESGetNormSchedule(snes->pc,&normschedule);

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:       SNESGetFunction(snes->pc,&FPC,NULL,NULL);
114:       if (FPC) {
115:         if (fnorm) {VecNorm(FPC,NORM_2,fnorm);}
116:         VecCopy(FPC,F);
117:       } else {
118:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function");
119:       }
120:     } else {
121:       SNESGetSolution(snes->pc,&XPC);
122:       if (XPC) {
123:         SNESComputeFunction(snes->pc,XPC,F);
124:         if (fnorm) {VecNorm(F,NORM_2,fnorm);}
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:   return(0);
129: }