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 /*@ 9*37703b41SBarry Smith SNESApplyPC - Calls SNESSolve() on preconditioner for the SNES 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 ierr = VecCopy(x,y);CHKERRQ(ierr); 4632f3f7c2SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr); 4732f3f7c2SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,y);CHKERRQ(ierr); 4832f3f7c2SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,x,y,0);CHKERRQ(ierr); 4932f3f7c2SPeter Brune ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr); 5032f3f7c2SPeter Brune PetscFunctionReturn(0); 5132f3f7c2SPeter Brune } 5232f3f7c2SPeter Brune PetscFunctionReturn(0); 5332f3f7c2SPeter Brune } 5432f3f7c2SPeter Brune 5532f3f7c2SPeter Brune #undef __FUNCT__ 5632f3f7c2SPeter Brune #define __FUNCT__ "SNESComputeFunctionDefaultPC" 57ed07d7d7SPeter Brune PetscErrorCode SNESComputeFunctionDefaultPC(SNES snes,Vec X,Vec F) { 5832f3f7c2SPeter Brune /* This is to be used as an argument to SNESMF -- NOT as a "function" */ 59b7281c8aSPeter Brune SNESConvergedReason reason; 6032f3f7c2SPeter Brune PetscErrorCode ierr; 61b7281c8aSPeter Brune 6232f3f7c2SPeter Brune PetscFunctionBegin; 63b7281c8aSPeter Brune if (snes->pc) { 646950c336SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 65b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 66b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 67b7281c8aSPeter Brune ierr = SNESSetFunctionDomainError(snes);CHKERRQ(ierr); 68b7281c8aSPeter Brune } 69b7281c8aSPeter Brune } else { 70b7281c8aSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 71b7281c8aSPeter Brune } 7232f3f7c2SPeter Brune PetscFunctionReturn(0); 7332f3f7c2SPeter Brune } 74ddd40ce5SPeter Brune 75ddd40ce5SPeter Brune #undef __FUNCT__ 76ddd40ce5SPeter Brune #define __FUNCT__ "SNESGetPCFunction" 77ddd40ce5SPeter Brune /*@ 78ddd40ce5SPeter Brune SNESGetPCFunction - Gets the function from a preconditioner after SNESSolve() has been called. 79ddd40ce5SPeter Brune 80ddd40ce5SPeter Brune Collective on SNES 81ddd40ce5SPeter Brune 82ddd40ce5SPeter Brune Input Parameters: 83ddd40ce5SPeter Brune . snes - the SNES context 84ddd40ce5SPeter Brune 85ddd40ce5SPeter Brune Output Parameter: 86ddd40ce5SPeter Brune . F - function vector 87ddd40ce5SPeter Brune . fnorm - the norm of F 88ddd40ce5SPeter Brune 89ddd40ce5SPeter Brune Level: developer 90ddd40ce5SPeter Brune 91ddd40ce5SPeter Brune .keywords: SNES, nonlinear, function 92ddd40ce5SPeter Brune 93ddd40ce5SPeter Brune .seealso: SNESGetPC(),SNESSetPC(),SNESComputeFunction(),SNESApplyPC(),SNESSolve() 94ddd40ce5SPeter Brune @*/ 95ddd40ce5SPeter Brune PetscErrorCode SNESGetPCFunction(SNES snes,Vec F,PetscReal *fnorm) 96ddd40ce5SPeter Brune { 97ddd40ce5SPeter Brune PetscErrorCode ierr; 98ddd40ce5SPeter Brune PCSide npcside; 99ddd40ce5SPeter Brune SNESFunctionType functype; 100ddd40ce5SPeter Brune SNESNormSchedule normschedule; 101ddd40ce5SPeter Brune Vec FPC,XPC; 102ddd40ce5SPeter Brune 103ddd40ce5SPeter Brune PetscFunctionBegin; 104ddd40ce5SPeter Brune if (snes->pc) { 105ddd40ce5SPeter Brune ierr = SNESGetPCSide(snes->pc,&npcside);CHKERRQ(ierr); 106ddd40ce5SPeter Brune ierr = SNESGetFunctionType(snes->pc,&functype);CHKERRQ(ierr); 107ddd40ce5SPeter Brune ierr = SNESGetNormSchedule(snes->pc,&normschedule);CHKERRQ(ierr); 108ddd40ce5SPeter Brune 109ddd40ce5SPeter Brune /* check if the function is valid based upon how the inner solver is preconditioned */ 110ddd40ce5SPeter Brune if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) { 111ddd40ce5SPeter Brune ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr); 1121d3ed18cSPeter Brune if (FPC) { 113ddd40ce5SPeter Brune if (fnorm) {ierr = SNESGetFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);} 114ddd40ce5SPeter Brune ierr = VecCopy(FPC,F);CHKERRQ(ierr); 115ddd40ce5SPeter Brune } else { 116ddd40ce5SPeter Brune SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function"); 117ddd40ce5SPeter Brune } 118ddd40ce5SPeter Brune } else { 119ddd40ce5SPeter Brune ierr = SNESGetSolution(snes->pc,&XPC);CHKERRQ(ierr); 120ddd40ce5SPeter Brune if (XPC) { 121ddd40ce5SPeter Brune ierr = SNESComputeFunction(snes->pc,XPC,F);CHKERRQ(ierr); 122ddd40ce5SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 123*37703b41SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution"); 124ddd40ce5SPeter Brune } 125*37703b41SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set"); 126ddd40ce5SPeter Brune PetscFunctionReturn(0); 127ddd40ce5SPeter Brune } 128