14800dd8cSBarry Smith 2c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 34800dd8cSBarry Smith 4971273eeSBarry Smith typedef struct { 5971273eeSBarry Smith SNES snes; 6df8705c3SBarry Smith /* Information on the regular SNES convergence test; which may have been user provided */ 7df8705c3SBarry Smith PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 8df8705c3SBarry Smith PetscErrorCode (*convdestroy)(void*); 9df8705c3SBarry Smith void *convctx; 10971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 11971273eeSBarry Smith 12470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 13df60cc22SBarry Smith { 14971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 15971273eeSBarry Smith SNES snes = ctx->snes; 1604d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 17df60cc22SBarry Smith Vec x; 18064f8208SBarry Smith PetscReal nrm; 19dfbe8321SBarry Smith PetscErrorCode ierr; 20df60cc22SBarry Smith 213a40ed3dSBarry Smith PetscFunctionBegin; 22df8705c3SBarry Smith ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr); 23329f5518SBarry Smith if (*reason) { 24df8705c3SBarry Smith ierr = PetscInfo2(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 259b04593eSLois Curfman McInnes } 26a935fc98SLois Curfman McInnes /* Determine norm of solution */ 2778b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 28064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 29064f8208SBarry Smith if (nrm >= neP->delta) { 3057622a8eSBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 31329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 32df60cc22SBarry Smith } 333a40ed3dSBarry Smith PetscFunctionReturn(0); 34df60cc22SBarry Smith } 3582bf6240SBarry Smith 36470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 37971273eeSBarry Smith { 38971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 39971273eeSBarry Smith PetscErrorCode ierr; 40971273eeSBarry Smith 41971273eeSBarry Smith PetscFunctionBegin; 42df8705c3SBarry Smith ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr); 43971273eeSBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 44971273eeSBarry Smith PetscFunctionReturn(0); 45971273eeSBarry Smith } 46971273eeSBarry Smith 4785385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 4885385478SLisandro Dalcin /* 49470880abSPatrick Sanan SNESTR_Converged_Private -test convergence JUST for 5085385478SLisandro Dalcin the trust region tolerance. 5185385478SLisandro Dalcin 5285385478SLisandro Dalcin */ 53470880abSPatrick Sanan static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 5485385478SLisandro Dalcin { 5504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 5685385478SLisandro Dalcin PetscErrorCode ierr; 5785385478SLisandro Dalcin 5885385478SLisandro Dalcin PetscFunctionBegin; 5985385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 6085385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 6157622a8eSBarry Smith ierr = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 621c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 63e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 6485385478SLisandro Dalcin ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 6585385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 6685385478SLisandro Dalcin } 6785385478SLisandro Dalcin PetscFunctionReturn(0); 6885385478SLisandro Dalcin } 6985385478SLisandro Dalcin 707cb011f5SBarry Smith /*@C 71*c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 72*c9368356SGlenn Hammond Allows the user a chance to change or override the decision of the line search routine. 73*c9368356SGlenn Hammond 74*c9368356SGlenn Hammond Logically Collective on snes 75*c9368356SGlenn Hammond 76*c9368356SGlenn Hammond Input Parameters: 77*c9368356SGlenn Hammond + snes - the nonlinear solver object 78*c9368356SGlenn Hammond . func - [optional] function evaluation routine, see SNESNewtonTRPreCheck() for the calling sequence 79*c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 80*c9368356SGlenn Hammond 81*c9368356SGlenn Hammond Level: intermediate 82*c9368356SGlenn Hammond 83*c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver. 84*c9368356SGlenn Hammond 85*c9368356SGlenn Hammond .seealso: SNESNewtonTRPreCheck(), SNESNewtonTRGetPreCheck() 86*c9368356SGlenn Hammond @*/ 87*c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx) 88*c9368356SGlenn Hammond { 89*c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 90*c9368356SGlenn Hammond 91*c9368356SGlenn Hammond PetscFunctionBegin; 92*c9368356SGlenn Hammond PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 93*c9368356SGlenn Hammond if (func) tr->precheck = func; 94*c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 95*c9368356SGlenn Hammond PetscFunctionReturn(0); 96*c9368356SGlenn Hammond } 97*c9368356SGlenn Hammond 98*c9368356SGlenn Hammond /*@C 99*c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 100*c9368356SGlenn Hammond 101*c9368356SGlenn Hammond Not collective 102*c9368356SGlenn Hammond 103*c9368356SGlenn Hammond Input Parameter: 104*c9368356SGlenn Hammond . snes - the nonlinear solver context 105*c9368356SGlenn Hammond 106*c9368356SGlenn Hammond Output Parameters: 107*c9368356SGlenn Hammond + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck() 108*c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 109*c9368356SGlenn Hammond 110*c9368356SGlenn Hammond Level: intermediate 111*c9368356SGlenn Hammond 112*c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRPreCheck() 113*c9368356SGlenn Hammond @*/ 114*c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx) 115*c9368356SGlenn Hammond { 116*c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 117*c9368356SGlenn Hammond 118*c9368356SGlenn Hammond PetscFunctionBegin; 119*c9368356SGlenn Hammond PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 120*c9368356SGlenn Hammond if (func) *func = tr->precheck; 121*c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 122*c9368356SGlenn Hammond PetscFunctionReturn(0); 123*c9368356SGlenn Hammond } 124*c9368356SGlenn Hammond 125*c9368356SGlenn Hammond /*@C 1267cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1277cb011f5SBarry Smith function evaluation. Allows the user a chance to change or override the decision of the line search routine 1287cb011f5SBarry Smith 1297cb011f5SBarry Smith Logically Collective on snes 1307cb011f5SBarry Smith 1317cb011f5SBarry Smith Input Parameters: 1327cb011f5SBarry Smith + snes - the nonlinear solver object 1337cb011f5SBarry Smith . func - [optional] function evaluation routine, see SNESNewtonTRPostCheck() for the calling sequence 1347cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1357cb011f5SBarry Smith 1367cb011f5SBarry Smith Level: intermediate 1377cb011f5SBarry Smith 138*c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in 1397cb011f5SBarry Smith SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 1407cb011f5SBarry Smith 1417cb011f5SBarry Smith .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck() 1427cb011f5SBarry Smith @*/ 143*c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 1447cb011f5SBarry Smith { 1457cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 1467cb011f5SBarry Smith 1477cb011f5SBarry Smith PetscFunctionBegin; 1487cb011f5SBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1497cb011f5SBarry Smith if (func) tr->postcheck = func; 1507cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1517cb011f5SBarry Smith PetscFunctionReturn(0); 1527cb011f5SBarry Smith } 1537cb011f5SBarry Smith 1547cb011f5SBarry Smith /*@C 1557cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1567cb011f5SBarry Smith 1577cb011f5SBarry Smith Not collective 1587cb011f5SBarry Smith 1597cb011f5SBarry Smith Input Parameter: 1607cb011f5SBarry Smith . snes - the nonlinear solver context 1617cb011f5SBarry Smith 1627cb011f5SBarry Smith Output Parameters: 1637cb011f5SBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck() 1647cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1657cb011f5SBarry Smith 1667cb011f5SBarry Smith Level: intermediate 1677cb011f5SBarry Smith 1687cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck() 1697cb011f5SBarry Smith @*/ 170*c9368356SGlenn Hammond PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 1717cb011f5SBarry Smith { 1727cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 1737cb011f5SBarry Smith 1747cb011f5SBarry Smith PetscFunctionBegin; 1757cb011f5SBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1767cb011f5SBarry Smith if (func) *func = tr->postcheck; 1777cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 1787cb011f5SBarry Smith PetscFunctionReturn(0); 1797cb011f5SBarry Smith } 1807cb011f5SBarry Smith 1817cb011f5SBarry Smith /*@C 182*c9368356SGlenn Hammond SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR 183*c9368356SGlenn Hammond 184*c9368356SGlenn Hammond Logically Collective on snes 185*c9368356SGlenn Hammond 186*c9368356SGlenn Hammond Input Parameters: 187*c9368356SGlenn Hammond + snes - the solver 188*c9368356SGlenn Hammond . X - The last solution 189*c9368356SGlenn Hammond - Y - The step direction 190*c9368356SGlenn Hammond 191*c9368356SGlenn Hammond Output Parameters: 192*c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 193*c9368356SGlenn Hammond 194*c9368356SGlenn Hammond Level: developer 195*c9368356SGlenn Hammond 196*c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck() 197*c9368356SGlenn Hammond @*/ 198*c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y) 199*c9368356SGlenn Hammond { 200*c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 201*c9368356SGlenn Hammond PetscErrorCode ierr; 202*c9368356SGlenn Hammond 203*c9368356SGlenn Hammond PetscFunctionBegin; 204*c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 205*c9368356SGlenn Hammond if (tr->precheck) { 206*c9368356SGlenn Hammond ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr); 207*c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes,*changed_Y,4); 208*c9368356SGlenn Hammond } 209*c9368356SGlenn Hammond PetscFunctionReturn(0); 210*c9368356SGlenn Hammond } 211*c9368356SGlenn Hammond 212*c9368356SGlenn Hammond /*@C 2137cb011f5SBarry Smith SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation 2147cb011f5SBarry Smith 2157cb011f5SBarry Smith Logically Collective on snes 2167cb011f5SBarry Smith 2177cb011f5SBarry Smith Input Parameters: 2187cb011f5SBarry Smith + snes - the solver. X - The last solution 2197cb011f5SBarry Smith . Y - The full step direction 2207cb011f5SBarry Smith - W - The updated solution, W = X + tao*Y for some tao 2217cb011f5SBarry Smith 2227cb011f5SBarry Smith Output Parameters: 2237cb011f5SBarry Smith . changed_W - Indicator if the new candidate solution W has been changed. 2247cb011f5SBarry Smith 2257cb011f5SBarry Smith Notes: 2267cb011f5SBarry Smith If Y is changed then W is recomputed as X + tao*Y where tao is the current update value in SNESNEWTONTR 2277cb011f5SBarry Smith 2287cb011f5SBarry Smith Level: developer 2297cb011f5SBarry Smith 2307cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck() 2317cb011f5SBarry Smith @*/ 232*c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 2337cb011f5SBarry Smith { 2347cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 2357cb011f5SBarry Smith PetscErrorCode ierr; 2367cb011f5SBarry Smith 2377cb011f5SBarry Smith PetscFunctionBegin; 238*c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2397cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2407cb011f5SBarry Smith if (tr->postcheck) { 241*c9368356SGlenn Hammond ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr); 242*c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes,*changed_Y,5); 2437cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes,*changed_W,6); 2447cb011f5SBarry Smith } 2457cb011f5SBarry Smith PetscFunctionReturn(0); 2467cb011f5SBarry Smith } 24785385478SLisandro Dalcin 248df60cc22SBarry Smith /* 24904d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 250ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 2514800dd8cSBarry Smith 252ddbbdb52SLois Curfman McInnes 2534800dd8cSBarry Smith */ 25404d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 2554800dd8cSBarry Smith { 25604d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 257*c9368356SGlenn 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; 266df8705c3SBarry Smith PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 2674800dd8cSBarry Smith 2683a40ed3dSBarry Smith PetscFunctionBegin; 2696c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 270c579b300SPatrick Farrell 271fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 272fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 27339e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 274fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 275fbe28522SBarry Smith G = snes->work[1]; 2766b5873e3SBarry Smith Ytmp = snes->work[2]; 277*c9368356SGlenn Hammond W = snes->work[3]; 2784800dd8cSBarry Smith 279e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 2804c49b128SBarry Smith snes->iter = 0; 281e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2824800dd8cSBarry Smith 283df8705c3SBarry Smith /* Set the linear stopping criteria to use the More' trick. */ 284df8705c3SBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 285df8705c3SBarry Smith ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr); 286df8705c3SBarry Smith if (convtest != SNESTR_KSPConverged_Private) { 287df8705c3SBarry Smith ierr = PetscNew(&ctx);CHKERRQ(ierr); 288df8705c3SBarry Smith ctx->snes = snes; 289df8705c3SBarry Smith ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 290df8705c3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr); 291df8705c3SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr); 292df8705c3SBarry Smith } 293df8705c3SBarry Smith 294e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2955334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 2961aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 297e4ed7901SPeter Brune 298cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 299422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 3008b84755bSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 301e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 302fbe28522SBarry Smith snes->norm = fnorm; 303e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 3048b84755bSBarry Smith delta = xnorm ? neP->delta0*xnorm : neP->delta0; 3054800dd8cSBarry Smith neP->delta = delta; 306a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3077a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 308b37302c6SLois Curfman McInnes 30985385478SLisandro Dalcin /* test convergence */ 31085385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 31185385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 3123f149594SLisandro Dalcin 313df60cc22SBarry Smith 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); 3271aa26658SKarl Rupp 328329f5518SBarry Smith snes->linear_its += lits; 3291aa26658SKarl Rupp 330ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 331064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 332064f8208SBarry Smith norm1 = nrm; 3336b5873e3SBarry Smith while (1) { 334*c9368356SGlenn 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; 34457622a8eSBarry Smith ierr = PetscInfo1(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 } 3531592aa04SStefano Zampini ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 354*c9368356SGlenn Hammond /* PreCheck() allows for updates to Y prior to W <- X - Y */ 355*c9368356SGlenn Hammond ierr = SNESNewtonTRPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr); 356*c9368356SGlenn Hammond ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); /* W <- X - Y */ 357*c9368356SGlenn Hammond ierr = SNESNewtonTRPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 358*c9368356SGlenn Hammond if (changed_y) ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); 359*c9368356SGlenn Hammond ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /* F(X) */ 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; 36957622a8eSBarry Smith ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 37057622a8eSBarry Smith ierr = PetscInfo3(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); 3756b5873e3SBarry Smith /* check to see if progress is hopeless */ 376ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 377470880abSPatrick Sanan ierr = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 37885385478SLisandro Dalcin if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);} 3797cb011f5SBarry Smith if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 380184914b5SBarry Smith if (reason) { 38152392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 3827a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 383a7cc72afSBarry Smith breakout = PETSC_TRUE; 384454a90a3SBarry Smith break; 38552392280SLois Curfman McInnes } 38650ffb88aSMatthew Knepley snes->numFailures++; 3874800dd8cSBarry Smith } 3881acabf8cSLois Curfman McInnes if (!breakout) { 38985385478SLisandro Dalcin /* Update function and solution vectors */ 3904800dd8cSBarry Smith fnorm = gnorm; 39185385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 392*c9368356SGlenn Hammond ierr = VecCopy(W,X);CHKERRQ(ierr); 39385385478SLisandro Dalcin /* Monitor convergence */ 394e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 395c509a162SBarry Smith snes->iter = i+1; 396fbe28522SBarry Smith snes->norm = fnorm; 397c1e67a49SFande Kong snes->xnorm = xnorm; 398c1e67a49SFande Kong snes->ynorm = ynorm; 399e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 400a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 4017a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 40285385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 403ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 404e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 405e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 4063f149594SLisandro Dalcin if (reason) break; 4071aa26658SKarl Rupp } else break; 40838442cffSBarry Smith } 40952392280SLois Curfman McInnes if (i == maxits) { 410ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 41185385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 41252392280SLois Curfman McInnes } 413e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 414184914b5SBarry Smith snes->reason = reason; 415e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 4163a40ed3dSBarry Smith PetscFunctionReturn(0); 4174800dd8cSBarry Smith } 4184800dd8cSBarry Smith /*------------------------------------------------------------*/ 41904d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 4204800dd8cSBarry Smith { 421dfbe8321SBarry Smith PetscErrorCode ierr; 4223a40ed3dSBarry Smith 4233a40ed3dSBarry Smith PetscFunctionBegin; 424*c9368356SGlenn Hammond ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 4256cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 4263a40ed3dSBarry Smith PetscFunctionReturn(0); 4274800dd8cSBarry Smith } 4286b8b9a38SLisandro Dalcin 42904d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 4306b8b9a38SLisandro Dalcin { 4316b8b9a38SLisandro Dalcin 4326b8b9a38SLisandro Dalcin PetscFunctionBegin; 4336b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 4346b8b9a38SLisandro Dalcin } 4356b8b9a38SLisandro Dalcin 43604d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 4374800dd8cSBarry Smith { 438dfbe8321SBarry Smith PetscErrorCode ierr; 4395baf8537SBarry Smith 4403a40ed3dSBarry Smith PetscFunctionBegin; 44104d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 442606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 4433a40ed3dSBarry Smith PetscFunctionReturn(0); 4444800dd8cSBarry Smith } 4454800dd8cSBarry Smith /*------------------------------------------------------------*/ 4464800dd8cSBarry Smith 4474416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes) 4484800dd8cSBarry Smith { 44904d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 450dfbe8321SBarry Smith PetscErrorCode ierr; 4514800dd8cSBarry Smith 4523a40ed3dSBarry Smith PetscFunctionBegin; 453e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 45494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 45594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 45694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 45794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 45894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 45994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 46094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 46194ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 462b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4633a40ed3dSBarry Smith PetscFunctionReturn(0); 4644800dd8cSBarry Smith } 4654800dd8cSBarry Smith 46604d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 467a935fc98SLois Curfman McInnes { 46804d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 469dfbe8321SBarry Smith PetscErrorCode ierr; 470ace3abfcSBarry Smith PetscBool iascii; 471a935fc98SLois Curfman McInnes 4723a40ed3dSBarry Smith PetscFunctionBegin; 473251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 47432077d6dSBarry Smith if (iascii) { 47554fe5c21SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr); 47657622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 47757622a8eSBarry 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); 47819bcc07fSBarry Smith } 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 480a935fc98SLois Curfman McInnes } 48152392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 482ebe3b25bSBarry Smith /*MC 48304d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 484ebe3b25bSBarry Smith 485ebe3b25bSBarry Smith Options Database: 4868e434772SBarry Smith + -snes_trtol <tol> - trust region tolerance 4878e434772SBarry Smith . -snes_tr_mu <mu> - trust region parameter 4888e434772SBarry Smith . -snes_tr_eta <eta> - trust region parameter 4898e434772SBarry Smith . -snes_tr_sigma <sigma> - trust region parameter 4908b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 4918e434772SBarry Smith . -snes_tr_delta1 <delta1> - trust region parameter 4928e434772SBarry Smith . -snes_tr_delta2 <delta2> - trust region parameter 4938e434772SBarry Smith - -snes_tr_delta3 <delta3> - trust region parameter 494ebe3b25bSBarry Smith 495ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 496ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 497ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 498ebe3b25bSBarry Smith 499ee3001cbSBarry Smith Level: intermediate 500ee3001cbSBarry Smith 50104d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 502ebe3b25bSBarry Smith 503ebe3b25bSBarry Smith M*/ 5048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 5054800dd8cSBarry Smith { 50604d7464bSBarry Smith SNES_NEWTONTR *neP; 507dfbe8321SBarry Smith PetscErrorCode ierr; 5084800dd8cSBarry Smith 5093a40ed3dSBarry Smith PetscFunctionBegin; 51004d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 51104d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 51204d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 51304d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 51404d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 51504d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 516fbe28522SBarry Smith 51742f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 518efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 51942f4f86dSBarry Smith 5204fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5214fc747eaSLawrence Mitchell 522b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 523fbe28522SBarry Smith snes->data = (void*)neP; 524fbe28522SBarry Smith neP->mu = 0.25; 525fbe28522SBarry Smith neP->eta = 0.75; 526fbe28522SBarry Smith neP->delta = 0.0; 527fbe28522SBarry Smith neP->delta0 = 0.2; 528fbe28522SBarry Smith neP->delta1 = 0.3; 529fbe28522SBarry Smith neP->delta2 = 0.75; 530fbe28522SBarry Smith neP->delta3 = 2.0; 531fbe28522SBarry Smith neP->sigma = 0.0001; 532ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 53375567043SBarry Smith neP->rnorm0 = 0.0; 53475567043SBarry Smith neP->ttol = 0.0; 5353a40ed3dSBarry Smith PetscFunctionReturn(0); 5364800dd8cSBarry Smith } 53782bf6240SBarry Smith 538