1c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 24800dd8cSBarry Smith 3971273eeSBarry Smith typedef struct { 4971273eeSBarry Smith SNES snes; 5df8705c3SBarry Smith /* Information on the regular SNES convergence test; which may have been user provided */ 6df8705c3SBarry Smith PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 7df8705c3SBarry Smith PetscErrorCode (*convdestroy)(void*); 8df8705c3SBarry Smith void *convctx; 9971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 10971273eeSBarry Smith 11470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 12df60cc22SBarry Smith { 13971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 14971273eeSBarry Smith SNES snes = ctx->snes; 1504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 16df60cc22SBarry Smith Vec x; 17064f8208SBarry Smith PetscReal nrm; 18dfbe8321SBarry Smith PetscErrorCode ierr; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 21df8705c3SBarry Smith ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr); 22329f5518SBarry Smith if (*reason) { 237d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 249b04593eSLois Curfman McInnes } 25a935fc98SLois Curfman McInnes /* Determine norm of solution */ 269e5d0892SLisandro Dalcin ierr = KSPBuildSolution(ksp,NULL,&x);CHKERRQ(ierr); 27064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 28064f8208SBarry Smith if (nrm >= neP->delta) { 297d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 30329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 31df60cc22SBarry Smith } 323a40ed3dSBarry Smith PetscFunctionReturn(0); 33df60cc22SBarry Smith } 3482bf6240SBarry Smith 35470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 36971273eeSBarry Smith { 37971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 38971273eeSBarry Smith PetscErrorCode ierr; 39971273eeSBarry Smith 40971273eeSBarry Smith PetscFunctionBegin; 41df8705c3SBarry Smith ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr); 42971273eeSBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 43971273eeSBarry Smith PetscFunctionReturn(0); 44971273eeSBarry Smith } 45971273eeSBarry Smith 4685385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 4785385478SLisandro Dalcin /* 48470880abSPatrick Sanan SNESTR_Converged_Private -test convergence JUST for 4985385478SLisandro Dalcin the trust region tolerance. 5085385478SLisandro Dalcin 5185385478SLisandro Dalcin */ 52470880abSPatrick Sanan static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 5385385478SLisandro Dalcin { 5404d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 5585385478SLisandro Dalcin PetscErrorCode ierr; 5685385478SLisandro Dalcin 5785385478SLisandro Dalcin PetscFunctionBegin; 5885385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 5985385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 607d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 611c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 62e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 637d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 6485385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 6585385478SLisandro Dalcin } 6685385478SLisandro Dalcin PetscFunctionReturn(0); 6785385478SLisandro Dalcin } 6885385478SLisandro Dalcin 697cb011f5SBarry Smith /*@C 70c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 71c9368356SGlenn Hammond Allows the user a chance to change or override the decision of the line search routine. 72c9368356SGlenn Hammond 73c9368356SGlenn Hammond Logically Collective on snes 74c9368356SGlenn Hammond 75c9368356SGlenn Hammond Input Parameters: 76c9368356SGlenn Hammond + snes - the nonlinear solver object 77c9368356SGlenn Hammond . func - [optional] function evaluation routine, see SNESNewtonTRPreCheck() for the calling sequence 78c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 79c9368356SGlenn Hammond 80c9368356SGlenn Hammond Level: intermediate 81c9368356SGlenn Hammond 82c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver. 83c9368356SGlenn Hammond 843312a946SBarry Smith .seealso: SNESNewtonTRPreCheck(), SNESNewtonTRGetPreCheck(), SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck() 85c9368356SGlenn Hammond @*/ 86c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx) 87c9368356SGlenn Hammond { 88c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 89c9368356SGlenn Hammond 90c9368356SGlenn Hammond PetscFunctionBegin; 91c9368356SGlenn Hammond PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 92c9368356SGlenn Hammond if (func) tr->precheck = func; 93c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 94c9368356SGlenn Hammond PetscFunctionReturn(0); 95c9368356SGlenn Hammond } 96c9368356SGlenn Hammond 97c9368356SGlenn Hammond /*@C 98c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 99c9368356SGlenn Hammond 100c9368356SGlenn Hammond Not collective 101c9368356SGlenn Hammond 102c9368356SGlenn Hammond Input Parameter: 103c9368356SGlenn Hammond . snes - the nonlinear solver context 104c9368356SGlenn Hammond 105c9368356SGlenn Hammond Output Parameters: 106c9368356SGlenn Hammond + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck() 107c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 108c9368356SGlenn Hammond 109c9368356SGlenn Hammond Level: intermediate 110c9368356SGlenn Hammond 111c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRPreCheck() 112c9368356SGlenn Hammond @*/ 113c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx) 114c9368356SGlenn Hammond { 115c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 116c9368356SGlenn Hammond 117c9368356SGlenn Hammond PetscFunctionBegin; 118c9368356SGlenn Hammond PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 119c9368356SGlenn Hammond if (func) *func = tr->precheck; 120c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 121c9368356SGlenn Hammond PetscFunctionReturn(0); 122c9368356SGlenn Hammond } 123c9368356SGlenn Hammond 124c9368356SGlenn Hammond /*@C 1257cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1267cb011f5SBarry Smith function evaluation. Allows the user a chance to change or override the decision of the line search routine 1277cb011f5SBarry Smith 1287cb011f5SBarry Smith Logically Collective on snes 1297cb011f5SBarry Smith 1307cb011f5SBarry Smith Input Parameters: 1317cb011f5SBarry Smith + snes - the nonlinear solver object 1327cb011f5SBarry Smith . func - [optional] function evaluation routine, see SNESNewtonTRPostCheck() for the calling sequence 1337cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1347cb011f5SBarry Smith 1357cb011f5SBarry Smith Level: intermediate 1367cb011f5SBarry Smith 137c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in 1387cb011f5SBarry Smith SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 1397cb011f5SBarry Smith 1407cb011f5SBarry Smith .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck() 1417cb011f5SBarry Smith @*/ 142c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 1437cb011f5SBarry Smith { 1447cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 1457cb011f5SBarry Smith 1467cb011f5SBarry Smith PetscFunctionBegin; 1477cb011f5SBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1487cb011f5SBarry Smith if (func) tr->postcheck = func; 1497cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1507cb011f5SBarry Smith PetscFunctionReturn(0); 1517cb011f5SBarry Smith } 1527cb011f5SBarry Smith 1537cb011f5SBarry Smith /*@C 1547cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1557cb011f5SBarry Smith 1567cb011f5SBarry Smith Not collective 1577cb011f5SBarry Smith 1587cb011f5SBarry Smith Input Parameter: 1597cb011f5SBarry Smith . snes - the nonlinear solver context 1607cb011f5SBarry Smith 1617cb011f5SBarry Smith Output Parameters: 1627cb011f5SBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck() 1637cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1647cb011f5SBarry Smith 1657cb011f5SBarry Smith Level: intermediate 1667cb011f5SBarry Smith 1677cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck() 1687cb011f5SBarry Smith @*/ 169c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes,PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 1707cb011f5SBarry Smith { 1717cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 1727cb011f5SBarry Smith 1737cb011f5SBarry Smith PetscFunctionBegin; 1747cb011f5SBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1757cb011f5SBarry Smith if (func) *func = tr->postcheck; 1767cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 1777cb011f5SBarry Smith PetscFunctionReturn(0); 1787cb011f5SBarry Smith } 1797cb011f5SBarry Smith 1807cb011f5SBarry Smith /*@C 181c9368356SGlenn Hammond SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR 182c9368356SGlenn Hammond 183c9368356SGlenn Hammond Logically Collective on snes 184c9368356SGlenn Hammond 185c9368356SGlenn Hammond Input Parameters: 186c9368356SGlenn Hammond + snes - the solver 187c9368356SGlenn Hammond . X - The last solution 188c9368356SGlenn Hammond - Y - The step direction 189c9368356SGlenn Hammond 190c9368356SGlenn Hammond Output Parameters: 191c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 192c9368356SGlenn Hammond 193c9368356SGlenn Hammond Level: developer 194c9368356SGlenn Hammond 195c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck() 196c9368356SGlenn Hammond @*/ 197c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y) 198c9368356SGlenn Hammond { 199c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 200c9368356SGlenn Hammond PetscErrorCode ierr; 201c9368356SGlenn Hammond 202c9368356SGlenn Hammond PetscFunctionBegin; 203c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 204c9368356SGlenn Hammond if (tr->precheck) { 205c9368356SGlenn Hammond ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr); 206c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes,*changed_Y,4); 207c9368356SGlenn Hammond } 208c9368356SGlenn Hammond PetscFunctionReturn(0); 209c9368356SGlenn Hammond } 210c9368356SGlenn Hammond 211c9368356SGlenn Hammond /*@C 2127cb011f5SBarry Smith SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation 2137cb011f5SBarry Smith 2147cb011f5SBarry Smith Logically Collective on snes 2157cb011f5SBarry Smith 2167cb011f5SBarry Smith Input Parameters: 2176b867d5aSJose E. Roman + snes - the solver 2186b867d5aSJose E. Roman . X - The last solution 2197cb011f5SBarry Smith . Y - The full step direction 2203312a946SBarry Smith - W - The updated solution, W = X - Y 2217cb011f5SBarry Smith 2227cb011f5SBarry Smith Output Parameters: 2233312a946SBarry Smith + changed_Y - indicator if step has been changed 2243312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 2257cb011f5SBarry Smith 2267cb011f5SBarry Smith Notes: 2273312a946SBarry Smith If Y is changed then W is recomputed as X - Y 2287cb011f5SBarry Smith 2297cb011f5SBarry Smith Level: developer 2307cb011f5SBarry Smith 2317cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck() 2327cb011f5SBarry Smith @*/ 233c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 2347cb011f5SBarry Smith { 2357cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 2367cb011f5SBarry Smith PetscErrorCode ierr; 2377cb011f5SBarry Smith 2387cb011f5SBarry Smith PetscFunctionBegin; 239c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2407cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2417cb011f5SBarry Smith if (tr->postcheck) { 242c9368356SGlenn Hammond ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr); 243c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes,*changed_Y,5); 2447cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes,*changed_W,6); 2457cb011f5SBarry Smith } 2467cb011f5SBarry Smith PetscFunctionReturn(0); 2477cb011f5SBarry Smith } 24885385478SLisandro Dalcin 249df60cc22SBarry Smith /* 25004d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 251ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 2524800dd8cSBarry Smith 2534800dd8cSBarry Smith */ 25404d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 2554800dd8cSBarry Smith { 25604d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 257c9368356SGlenn Hammond Vec X,F,Y,G,Ytmp,W; 2586849ba73SBarry Smith PetscErrorCode ierr; 259a7cc72afSBarry Smith PetscInt maxits,i,lits; 26085385478SLisandro Dalcin PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 26179f36460SBarry Smith PetscScalar cnorm; 262df60cc22SBarry Smith KSP ksp; 2633f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 264df8705c3SBarry Smith PetscBool breakout = PETSC_FALSE; 265df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 2665e28bcb6Sprj- PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),(*convdestroy)(void*); 2675e28bcb6Sprj- void *convctx; 2684800dd8cSBarry Smith 2693a40ed3dSBarry Smith PetscFunctionBegin; 270*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 271c579b300SPatrick Farrell 272fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 273fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 27439e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 275fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 276fbe28522SBarry Smith G = snes->work[1]; 2776b5873e3SBarry Smith Ytmp = snes->work[2]; 278c9368356SGlenn Hammond W = snes->work[3]; 2794800dd8cSBarry Smith 280e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 2814c49b128SBarry Smith snes->iter = 0; 282e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2834800dd8cSBarry Smith 284df8705c3SBarry Smith /* Set the linear stopping criteria to use the More' trick. */ 285df8705c3SBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 2865e28bcb6Sprj- ierr = KSPGetConvergenceTest(ksp,&convtest,&convctx,&convdestroy);CHKERRQ(ierr); 287df8705c3SBarry Smith if (convtest != SNESTR_KSPConverged_Private) { 288df8705c3SBarry Smith ierr = PetscNew(&ctx);CHKERRQ(ierr); 289df8705c3SBarry Smith ctx->snes = snes; 290df8705c3SBarry Smith ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 291df8705c3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr); 292df8705c3SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr); 293df8705c3SBarry Smith } 294df8705c3SBarry Smith 295e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2965334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 2971aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 298e4ed7901SPeter Brune 299cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 300422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 301284fb49fSHeeho Park ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- || X || */ 302e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 303fbe28522SBarry Smith snes->norm = fnorm; 304e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 3058b84755bSBarry Smith delta = xnorm ? neP->delta0*xnorm : neP->delta0; 3064800dd8cSBarry Smith neP->delta = delta; 307a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3087a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 309b37302c6SLois Curfman McInnes 31085385478SLisandro Dalcin /* test convergence */ 31185385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 31285385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 3133f149594SLisandro Dalcin 3144800dd8cSBarry Smith for (i=0; i<maxits; i++) { 31599a96b7cSMatthew Knepley 31699a96b7cSMatthew Knepley /* Call general purpose update function */ 317e7788613SBarry Smith if (snes->ops->update) { 318e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 31999a96b7cSMatthew Knepley } 32099a96b7cSMatthew Knepley 32185385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 322d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 32307b62357SFande Kong SNESCheckJacobianDomainerror(snes); 32423ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 325d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 3263f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 327329f5518SBarry Smith snes->linear_its += lits; 3281aa26658SKarl Rupp 3297d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n",snes->iter,lits);CHKERRQ(ierr); 330064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 331064f8208SBarry Smith norm1 = nrm; 332284fb49fSHeeho Park 3336b5873e3SBarry Smith while (1) { 334c9368356SGlenn Hammond PetscBool changed_y; 3357cb011f5SBarry Smith PetscBool changed_w; 336393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 337064f8208SBarry Smith nrm = norm1; 3386b5873e3SBarry Smith 3392deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 340064f8208SBarry Smith if (nrm >= delta) { 341064f8208SBarry Smith nrm = delta/nrm; 342064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 343064f8208SBarry Smith cnorm = nrm; 3447d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 3452dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 346064f8208SBarry Smith nrm = gpnorm; 347fbe28522SBarry Smith ynorm = delta; 348fbe28522SBarry Smith } else { 349fbe28522SBarry Smith gpnorm = 0.0; 350ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 351064f8208SBarry Smith ynorm = nrm; 352fbe28522SBarry Smith } 353c9368356SGlenn Hammond /* PreCheck() allows for updates to Y prior to W <- X - Y */ 354c9368356SGlenn Hammond ierr = SNESNewtonTRPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr); 355c9368356SGlenn Hammond ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); /* W <- X - Y */ 356c9368356SGlenn Hammond ierr = SNESNewtonTRPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 35760d4fc61SSatish Balay if (changed_y) {ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);} 358c1cd17f6SGlenn Hammond ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 359284fb49fSHeeho Park ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /* F(X-Y) = G */ 360cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 3617551ef16SBarry Smith SNESCheckFunctionNorm(snes,gnorm); 3624800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 363fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 3644800dd8cSBarry Smith 3654800dd8cSBarry Smith /* Update size of trust region */ 3664800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 3674800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 3684800dd8cSBarry Smith else delta *= neP->delta3; 3697d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 3707d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 3711aa26658SKarl Rupp 3724800dd8cSBarry Smith neP->delta = delta; 3734800dd8cSBarry Smith if (rho > neP->sigma) break; 374ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 375284fb49fSHeeho Park 3766b5873e3SBarry Smith /* check to see if progress is hopeless */ 377ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 378470880abSPatrick Sanan ierr = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 37985385478SLisandro Dalcin if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);} 3807cb011f5SBarry Smith if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 381184914b5SBarry Smith if (reason) { 38252392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 3837a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 384a7cc72afSBarry Smith breakout = PETSC_TRUE; 385454a90a3SBarry Smith break; 38652392280SLois Curfman McInnes } 38750ffb88aSMatthew Knepley snes->numFailures++; 3884800dd8cSBarry Smith } 3891acabf8cSLois Curfman McInnes if (!breakout) { 39085385478SLisandro Dalcin /* Update function and solution vectors */ 3914800dd8cSBarry Smith fnorm = gnorm; 39285385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 393c9368356SGlenn Hammond ierr = VecCopy(W,X);CHKERRQ(ierr); 39485385478SLisandro Dalcin /* Monitor convergence */ 395e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 396c509a162SBarry Smith snes->iter = i+1; 397fbe28522SBarry Smith snes->norm = fnorm; 398c1e67a49SFande Kong snes->xnorm = xnorm; 399c1e67a49SFande Kong snes->ynorm = ynorm; 400e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 401a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 4027a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 40385385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 404ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 405e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) {ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);} 406e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 4073f149594SLisandro Dalcin if (reason) break; 4081aa26658SKarl Rupp } else break; 40938442cffSBarry Smith } 410284fb49fSHeeho Park 41152392280SLois Curfman McInnes if (i == maxits) { 4127d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",maxits);CHKERRQ(ierr); 41385385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 41452392280SLois Curfman McInnes } 415e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 416184914b5SBarry Smith snes->reason = reason; 417e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 4185e28bcb6Sprj- if (convtest != SNESTR_KSPConverged_Private) { 4195e28bcb6Sprj- ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 4205e28bcb6Sprj- ierr = PetscFree(ctx);CHKERRQ(ierr); 4215e28bcb6Sprj- ierr = KSPSetConvergenceTest(ksp,convtest,convctx,convdestroy);CHKERRQ(ierr); 4225e28bcb6Sprj- } 4233a40ed3dSBarry Smith PetscFunctionReturn(0); 4244800dd8cSBarry Smith } 425284fb49fSHeeho Park 4264800dd8cSBarry Smith /*------------------------------------------------------------*/ 42704d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 4284800dd8cSBarry Smith { 429dfbe8321SBarry Smith PetscErrorCode ierr; 4303a40ed3dSBarry Smith 4313a40ed3dSBarry Smith PetscFunctionBegin; 432c9368356SGlenn Hammond ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 4336cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 4343a40ed3dSBarry Smith PetscFunctionReturn(0); 4354800dd8cSBarry Smith } 4366b8b9a38SLisandro Dalcin 43704d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 4386b8b9a38SLisandro Dalcin { 4396b8b9a38SLisandro Dalcin 4406b8b9a38SLisandro Dalcin PetscFunctionBegin; 4416b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 4426b8b9a38SLisandro Dalcin } 4436b8b9a38SLisandro Dalcin 44404d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 4454800dd8cSBarry Smith { 446dfbe8321SBarry Smith PetscErrorCode ierr; 4475baf8537SBarry Smith 4483a40ed3dSBarry Smith PetscFunctionBegin; 44904d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 450606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 4513a40ed3dSBarry Smith PetscFunctionReturn(0); 4524800dd8cSBarry Smith } 4534800dd8cSBarry Smith /*------------------------------------------------------------*/ 4544800dd8cSBarry Smith 4554416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes) 4564800dd8cSBarry Smith { 45704d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 458dfbe8321SBarry Smith PetscErrorCode ierr; 4594800dd8cSBarry Smith 4603a40ed3dSBarry Smith PetscFunctionBegin; 461e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 46294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 46394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 46494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 46594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 46694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 46794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 46894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 46994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 470b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4713a40ed3dSBarry Smith PetscFunctionReturn(0); 4724800dd8cSBarry Smith } 4734800dd8cSBarry Smith 47404d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 475a935fc98SLois Curfman McInnes { 47604d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 477dfbe8321SBarry Smith PetscErrorCode ierr; 478ace3abfcSBarry Smith PetscBool iascii; 479a935fc98SLois Curfman McInnes 4803a40ed3dSBarry Smith PetscFunctionBegin; 481251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 48232077d6dSBarry Smith if (iascii) { 48354fe5c21SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr); 48457622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 48557622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr); 48619bcc07fSBarry Smith } 4873a40ed3dSBarry Smith PetscFunctionReturn(0); 488a935fc98SLois Curfman McInnes } 48952392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 490ebe3b25bSBarry Smith /*MC 49104d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 492ebe3b25bSBarry Smith 493ebe3b25bSBarry Smith Options Database: 4948e434772SBarry Smith + -snes_trtol <tol> - trust region tolerance 4958e434772SBarry Smith . -snes_tr_mu <mu> - trust region parameter 4968e434772SBarry Smith . -snes_tr_eta <eta> - trust region parameter 4978e434772SBarry Smith . -snes_tr_sigma <sigma> - trust region parameter 4988b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 4998e434772SBarry Smith . -snes_tr_delta1 <delta1> - trust region parameter 5008e434772SBarry Smith . -snes_tr_delta2 <delta2> - trust region parameter 5018e434772SBarry Smith - -snes_tr_delta3 <delta3> - trust region parameter 502ebe3b25bSBarry Smith 503ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 504ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 505ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 506ebe3b25bSBarry Smith 507ee3001cbSBarry Smith Level: intermediate 508ee3001cbSBarry Smith 50904d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 510ebe3b25bSBarry Smith 511ebe3b25bSBarry Smith M*/ 5128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 5134800dd8cSBarry Smith { 51404d7464bSBarry Smith SNES_NEWTONTR *neP; 515dfbe8321SBarry Smith PetscErrorCode ierr; 5164800dd8cSBarry Smith 5173a40ed3dSBarry Smith PetscFunctionBegin; 51804d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 51904d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 52004d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 52104d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 52204d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 52304d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 524fbe28522SBarry Smith 52542f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 526efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 52742f4f86dSBarry Smith 5284fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5294fc747eaSLawrence Mitchell 530b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 531fbe28522SBarry Smith snes->data = (void*)neP; 532fbe28522SBarry Smith neP->mu = 0.25; 533fbe28522SBarry Smith neP->eta = 0.75; 534fbe28522SBarry Smith neP->delta = 0.0; 535fbe28522SBarry Smith neP->delta0 = 0.2; 536fbe28522SBarry Smith neP->delta1 = 0.3; 537fbe28522SBarry Smith neP->delta2 = 0.75; 538fbe28522SBarry Smith neP->delta3 = 2.0; 539fbe28522SBarry Smith neP->sigma = 0.0001; 540ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 54175567043SBarry Smith neP->rnorm0 = 0.0; 54275567043SBarry Smith neP->ttol = 0.0; 5433a40ed3dSBarry Smith PetscFunctionReturn(0); 5444800dd8cSBarry Smith } 545