14800dd8cSBarry Smith 2c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 34800dd8cSBarry Smith 4971273eeSBarry Smith typedef struct { 5971273eeSBarry Smith void *ctx; 6971273eeSBarry Smith SNES snes; 7971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 8971273eeSBarry Smith 9fbe28522SBarry Smith /* 10df60cc22SBarry Smith This convergence test determines if the two norm of the 11df60cc22SBarry Smith solution lies outside the trust region, if so it halts. 12df60cc22SBarry Smith */ 134a2ae208SSatish Balay #undef __FUNCT__ 144b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private" 15971273eeSBarry Smith PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 16df60cc22SBarry Smith { 17971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 18971273eeSBarry Smith SNES snes = ctx->snes; 1904d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 20df60cc22SBarry Smith Vec x; 21064f8208SBarry Smith PetscReal nrm; 22dfbe8321SBarry Smith PetscErrorCode ierr; 23df60cc22SBarry Smith 243a40ed3dSBarry Smith PetscFunctionBegin; 258de4e6dcSJed Brown ierr = KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr); 26329f5518SBarry Smith if (*reason) { 2757622a8eSBarry Smith ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 289b04593eSLois Curfman McInnes } 29a935fc98SLois Curfman McInnes /* Determine norm of solution */ 3078b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 31064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 32064f8208SBarry Smith if (nrm >= neP->delta) { 3357622a8eSBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 34329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 35df60cc22SBarry Smith } 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37df60cc22SBarry Smith } 3882bf6240SBarry Smith 39971273eeSBarry Smith #undef __FUNCT__ 40971273eeSBarry Smith #define __FUNCT__ "SNES_TR_KSPConverged_Destroy" 41971273eeSBarry Smith PetscErrorCode SNES_TR_KSPConverged_Destroy(void *cctx) 42971273eeSBarry Smith { 43971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 44971273eeSBarry Smith PetscErrorCode ierr; 45971273eeSBarry Smith 46971273eeSBarry Smith PetscFunctionBegin; 478de4e6dcSJed Brown ierr = KSPConvergedDefaultDestroy(ctx->ctx);CHKERRQ(ierr); 48971273eeSBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 49971273eeSBarry Smith PetscFunctionReturn(0); 50971273eeSBarry Smith } 51971273eeSBarry Smith 5285385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 5385385478SLisandro Dalcin #undef __FUNCT__ 5485385478SLisandro Dalcin #define __FUNCT__ "SNES_TR_Converged_Private" 5585385478SLisandro Dalcin /* 5685385478SLisandro Dalcin SNES_TR_Converged_Private -test convergence JUST for 5785385478SLisandro Dalcin the trust region tolerance. 5885385478SLisandro Dalcin 5985385478SLisandro Dalcin */ 6085385478SLisandro Dalcin static PetscErrorCode SNES_TR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 6185385478SLisandro Dalcin { 6204d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 6385385478SLisandro Dalcin PetscErrorCode ierr; 6485385478SLisandro Dalcin 6585385478SLisandro Dalcin PetscFunctionBegin; 6685385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 6785385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 6857622a8eSBarry Smith ierr = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 6985385478SLisandro Dalcin *reason = SNES_CONVERGED_TR_DELTA; 7085385478SLisandro Dalcin } else if (snes->nfuncs >= snes->max_funcs) { 7185385478SLisandro Dalcin ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 7285385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 7385385478SLisandro Dalcin } 7485385478SLisandro Dalcin PetscFunctionReturn(0); 7585385478SLisandro Dalcin } 7685385478SLisandro Dalcin 7785385478SLisandro Dalcin 78df60cc22SBarry Smith /* 7904d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 80ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 814800dd8cSBarry Smith 82ddbbdb52SLois Curfman McInnes 834800dd8cSBarry Smith */ 844a2ae208SSatish Balay #undef __FUNCT__ 8504d7464bSBarry Smith #define __FUNCT__ "SNESSolve_NEWTONTR" 8604d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 874800dd8cSBarry Smith { 8804d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 8985385478SLisandro Dalcin Vec X,F,Y,G,Ytmp; 906849ba73SBarry Smith PetscErrorCode ierr; 91a7cc72afSBarry Smith PetscInt maxits,i,lits; 9285385478SLisandro Dalcin PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 9379f36460SBarry Smith PetscScalar cnorm; 94df60cc22SBarry Smith KSP ksp; 953f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 96ace3abfcSBarry Smith PetscBool conv = PETSC_FALSE,breakout = PETSC_FALSE; 974800dd8cSBarry Smith 983a40ed3dSBarry Smith PetscFunctionBegin; 99c579b300SPatrick Farrell 100c579b300SPatrick Farrell if (snes->xl || snes->xu || snes->ops->computevariablebounds) { 101c579b300SPatrick Farrell SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 102c579b300SPatrick Farrell } 103c579b300SPatrick Farrell 104fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 105fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 10639e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 107fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 108fbe28522SBarry Smith G = snes->work[1]; 1096b5873e3SBarry Smith Ytmp = snes->work[2]; 1104800dd8cSBarry Smith 111e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1124c49b128SBarry Smith snes->iter = 0; 113e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1144800dd8cSBarry Smith 115e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1165334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 1171aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 118e4ed7901SPeter Brune 119cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 120422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 121*8b84755bSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 122e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 123fbe28522SBarry Smith snes->norm = fnorm; 124e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 125*8b84755bSBarry Smith delta = xnorm ? neP->delta0*xnorm : neP->delta0; 1264800dd8cSBarry Smith neP->delta = delta; 127a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1287a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 129b37302c6SLois Curfman McInnes 13085385478SLisandro Dalcin /* test convergence */ 13185385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 13285385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 1333f149594SLisandro Dalcin 134a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 1350298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 1367c4af3dcSLois Curfman McInnes if (!conv) { 137971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 1383f149594SLisandro Dalcin ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 139b00a9115SJed Brown ierr = PetscNew(&ctx);CHKERRQ(ierr); 140971273eeSBarry Smith ctx->snes = snes; 1418de4e6dcSJed Brown ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 142971273eeSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr); 143ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 1447c4af3dcSLois Curfman McInnes } 145df60cc22SBarry Smith 1464800dd8cSBarry Smith for (i=0; i<maxits; i++) { 14799a96b7cSMatthew Knepley 14899a96b7cSMatthew Knepley /* Call general purpose update function */ 149e7788613SBarry Smith if (snes->ops->update) { 150e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 15199a96b7cSMatthew Knepley } 15299a96b7cSMatthew Knepley 15385385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 154d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 15523ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 156d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 1573f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1581aa26658SKarl Rupp 159329f5518SBarry Smith snes->linear_its += lits; 1601aa26658SKarl Rupp 161ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 162064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 163064f8208SBarry Smith norm1 = nrm; 1646b5873e3SBarry Smith while (1) { 165393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 166064f8208SBarry Smith nrm = norm1; 1676b5873e3SBarry Smith 1682deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 169064f8208SBarry Smith if (nrm >= delta) { 170064f8208SBarry Smith nrm = delta/nrm; 171064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 172064f8208SBarry Smith cnorm = nrm; 17357622a8eSBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 1742dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 175064f8208SBarry Smith nrm = gpnorm; 176fbe28522SBarry Smith ynorm = delta; 177fbe28522SBarry Smith } else { 178fbe28522SBarry Smith gpnorm = 0.0; 179ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 180064f8208SBarry Smith ynorm = nrm; 181fbe28522SBarry Smith } 18279f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 18385385478SLisandro Dalcin ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 1845334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 185cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1864800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 187fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1884800dd8cSBarry Smith 1894800dd8cSBarry Smith /* Update size of trust region */ 1904800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1914800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1924800dd8cSBarry Smith else delta *= neP->delta3; 19357622a8eSBarry Smith ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 19457622a8eSBarry Smith ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 1951aa26658SKarl Rupp 1964800dd8cSBarry Smith neP->delta = delta; 1974800dd8cSBarry Smith if (rho > neP->sigma) break; 198ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 1996b5873e3SBarry Smith /* check to see if progress is hopeless */ 200ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 20185385478SLisandro Dalcin ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 20285385478SLisandro Dalcin if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 203184914b5SBarry Smith if (reason) { 20452392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 2057a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 206a7cc72afSBarry Smith breakout = PETSC_TRUE; 207454a90a3SBarry Smith break; 20852392280SLois Curfman McInnes } 20950ffb88aSMatthew Knepley snes->numFailures++; 2104800dd8cSBarry Smith } 2111acabf8cSLois Curfman McInnes if (!breakout) { 21285385478SLisandro Dalcin /* Update function and solution vectors */ 2134800dd8cSBarry Smith fnorm = gnorm; 21485385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 21585385478SLisandro Dalcin ierr = VecCopy(Y,X);CHKERRQ(ierr); 21685385478SLisandro Dalcin /* Monitor convergence */ 217e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 218c509a162SBarry Smith snes->iter = i+1; 219fbe28522SBarry Smith snes->norm = fnorm; 220e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 221a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2227a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 22385385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 224ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 225e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 226e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 2273f149594SLisandro Dalcin if (reason) break; 2281aa26658SKarl Rupp } else break; 22938442cffSBarry Smith } 23052392280SLois Curfman McInnes if (i == maxits) { 231ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 23285385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 23352392280SLois Curfman McInnes } 234e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 235184914b5SBarry Smith snes->reason = reason; 236e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2373a40ed3dSBarry Smith PetscFunctionReturn(0); 2384800dd8cSBarry Smith } 2394800dd8cSBarry Smith /*------------------------------------------------------------*/ 2404a2ae208SSatish Balay #undef __FUNCT__ 24104d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONTR" 24204d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 2434800dd8cSBarry Smith { 244dfbe8321SBarry Smith PetscErrorCode ierr; 2453a40ed3dSBarry Smith 2463a40ed3dSBarry Smith PetscFunctionBegin; 247fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 2486cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2493a40ed3dSBarry Smith PetscFunctionReturn(0); 2504800dd8cSBarry Smith } 2516b8b9a38SLisandro Dalcin 2526b8b9a38SLisandro Dalcin #undef __FUNCT__ 25304d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONTR" 25404d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 2556b8b9a38SLisandro Dalcin { 2566b8b9a38SLisandro Dalcin 2576b8b9a38SLisandro Dalcin PetscFunctionBegin; 2586b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2596b8b9a38SLisandro Dalcin } 2606b8b9a38SLisandro Dalcin 2614a2ae208SSatish Balay #undef __FUNCT__ 26204d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONTR" 26304d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 2644800dd8cSBarry Smith { 265dfbe8321SBarry Smith PetscErrorCode ierr; 2665baf8537SBarry Smith 2673a40ed3dSBarry Smith PetscFunctionBegin; 26804d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 269606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2703a40ed3dSBarry Smith PetscFunctionReturn(0); 2714800dd8cSBarry Smith } 2724800dd8cSBarry Smith /*------------------------------------------------------------*/ 2734800dd8cSBarry Smith 2744a2ae208SSatish Balay #undef __FUNCT__ 27504d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONTR" 2768c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptions *PetscOptionsObject,SNES snes) 2774800dd8cSBarry Smith { 27804d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 279dfbe8321SBarry Smith PetscErrorCode ierr; 2804800dd8cSBarry Smith 2813a40ed3dSBarry Smith PetscFunctionBegin; 282e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 28394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 28494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 28594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 28694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 28794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 28894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 28994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 29094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 291b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2923a40ed3dSBarry Smith PetscFunctionReturn(0); 2934800dd8cSBarry Smith } 2944800dd8cSBarry Smith 2954a2ae208SSatish Balay #undef __FUNCT__ 29604d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONTR" 29704d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 298a935fc98SLois Curfman McInnes { 29904d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 300dfbe8321SBarry Smith PetscErrorCode ierr; 301ace3abfcSBarry Smith PetscBool iascii; 302a935fc98SLois Curfman McInnes 3033a40ed3dSBarry Smith PetscFunctionBegin; 304251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 30532077d6dSBarry Smith if (iascii) { 30657622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 30757622a8eSBarry 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); 30819bcc07fSBarry Smith } 3093a40ed3dSBarry Smith PetscFunctionReturn(0); 310a935fc98SLois Curfman McInnes } 31152392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 312ebe3b25bSBarry Smith /*MC 31304d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 314ebe3b25bSBarry Smith 315ebe3b25bSBarry Smith Options Database: 316ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 317ebe3b25bSBarry Smith . -snes_tr_mu <mu> 318ebe3b25bSBarry Smith . -snes_tr_eta <eta> 319ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 320*8b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 321ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 322ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 323ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 324ebe3b25bSBarry Smith 325ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 326ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 327ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 328ebe3b25bSBarry Smith 329ee3001cbSBarry Smith Level: intermediate 330ee3001cbSBarry Smith 33104d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 332ebe3b25bSBarry Smith 333ebe3b25bSBarry Smith M*/ 3344a2ae208SSatish Balay #undef __FUNCT__ 33504d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONTR" 3368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 3374800dd8cSBarry Smith { 33804d7464bSBarry Smith SNES_NEWTONTR *neP; 339dfbe8321SBarry Smith PetscErrorCode ierr; 3404800dd8cSBarry Smith 3413a40ed3dSBarry Smith PetscFunctionBegin; 34204d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 34304d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 34404d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 34504d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 34604d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 34704d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 348fbe28522SBarry Smith 34942f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 35042f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 35142f4f86dSBarry Smith 352b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 353fbe28522SBarry Smith snes->data = (void*)neP; 354fbe28522SBarry Smith neP->mu = 0.25; 355fbe28522SBarry Smith neP->eta = 0.75; 356fbe28522SBarry Smith neP->delta = 0.0; 357fbe28522SBarry Smith neP->delta0 = 0.2; 358fbe28522SBarry Smith neP->delta1 = 0.3; 359fbe28522SBarry Smith neP->delta2 = 0.75; 360fbe28522SBarry Smith neP->delta3 = 2.0; 361fbe28522SBarry Smith neP->sigma = 0.0001; 362ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 36375567043SBarry Smith neP->rnorm0 = 0.0; 36475567043SBarry Smith neP->ttol = 0.0; 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 3664800dd8cSBarry Smith } 36782bf6240SBarry Smith 368