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 the function that has been set with SNESSetFunction(). 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 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "SNESComputeFunctionDefaultPC" 59 PetscErrorCode SNESComputeFunctionDefaultPC(SNES snes,Vec X,Vec F) { 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 = SNESApplyPC(snes,X,NULL,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__ "SNESGetPCFunction" 79 /*@ 80 SNESGetPCFunction - 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: SNESGetPC(),SNESSetPC(),SNESComputeFunction(),SNESApplyPC(),SNESSolve() 96 @*/ 97 PetscErrorCode SNESGetPCFunction(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 = SNESGetPCSide(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 = SNESGetFunctionNorm(snes->pc,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 { 126 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution"); 127 } 128 } 129 } else { 130 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set"); 131 } 132 133 PetscFunctionReturn(0); 134 } 135