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; 97e4ed7901SPeter Brune PetscBool domainerror; 984800dd8cSBarry Smith 993a40ed3dSBarry Smith PetscFunctionBegin; 100*c579b300SPatrick Farrell 101*c579b300SPatrick Farrell if (snes->xl || snes->xu || snes->ops->computevariablebounds) { 102*c579b300SPatrick Farrell SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 103*c579b300SPatrick Farrell } 104*c579b300SPatrick Farrell 105fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 106fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 10739e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 108fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 109fbe28522SBarry Smith G = snes->work[1]; 1106b5873e3SBarry Smith Ytmp = snes->work[2]; 1114800dd8cSBarry Smith 112e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1134c49b128SBarry Smith snes->iter = 0; 114e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1154800dd8cSBarry Smith 116e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1175334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 118e4ed7901SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 119e4ed7901SPeter Brune if (domainerror) { 120e4ed7901SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 121e4ed7901SPeter Brune PetscFunctionReturn(0); 122e4ed7901SPeter Brune } 1231aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 124e4ed7901SPeter Brune 125cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 126189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 127189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 128189a9710SBarry Smith PetscFunctionReturn(0); 129189a9710SBarry Smith } 130e4ed7901SPeter Brune 131e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 132fbe28522SBarry Smith snes->norm = fnorm; 133e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1344800dd8cSBarry Smith delta = neP->delta0*fnorm; 1354800dd8cSBarry Smith neP->delta = delta; 136a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1377a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 138b37302c6SLois Curfman McInnes 13985385478SLisandro Dalcin /* test convergence */ 14085385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 14185385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 1423f149594SLisandro Dalcin 143a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 1440298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 1457c4af3dcSLois Curfman McInnes if (!conv) { 146971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 1473f149594SLisandro Dalcin ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 148b00a9115SJed Brown ierr = PetscNew(&ctx);CHKERRQ(ierr); 149971273eeSBarry Smith ctx->snes = snes; 1508de4e6dcSJed Brown ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 151971273eeSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr); 152ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 1537c4af3dcSLois Curfman McInnes } 154df60cc22SBarry Smith 1554800dd8cSBarry Smith for (i=0; i<maxits; i++) { 15699a96b7cSMatthew Knepley 15799a96b7cSMatthew Knepley /* Call general purpose update function */ 158e7788613SBarry Smith if (snes->ops->update) { 159e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 16099a96b7cSMatthew Knepley } 16199a96b7cSMatthew Knepley 16285385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 163d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 16423ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 165d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 1663f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1671aa26658SKarl Rupp 168329f5518SBarry Smith snes->linear_its += lits; 1691aa26658SKarl Rupp 170ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 171064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 172064f8208SBarry Smith norm1 = nrm; 1736b5873e3SBarry Smith while (1) { 174393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 175064f8208SBarry Smith nrm = norm1; 1766b5873e3SBarry Smith 1772deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 178064f8208SBarry Smith if (nrm >= delta) { 179064f8208SBarry Smith nrm = delta/nrm; 180064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 181064f8208SBarry Smith cnorm = nrm; 18257622a8eSBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 1832dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 184064f8208SBarry Smith nrm = gpnorm; 185fbe28522SBarry Smith ynorm = delta; 186fbe28522SBarry Smith } else { 187fbe28522SBarry Smith gpnorm = 0.0; 188ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 189064f8208SBarry Smith ynorm = nrm; 190fbe28522SBarry Smith } 19179f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 19285385478SLisandro Dalcin ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 1935334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 194cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1954800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 196fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1974800dd8cSBarry Smith 1984800dd8cSBarry Smith /* Update size of trust region */ 1994800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 2004800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 2014800dd8cSBarry Smith else delta *= neP->delta3; 20257622a8eSBarry Smith ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 20357622a8eSBarry Smith ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 2041aa26658SKarl Rupp 2054800dd8cSBarry Smith neP->delta = delta; 2064800dd8cSBarry Smith if (rho > neP->sigma) break; 207ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 2086b5873e3SBarry Smith /* check to see if progress is hopeless */ 209ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 21085385478SLisandro Dalcin ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 21185385478SLisandro Dalcin if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 212184914b5SBarry Smith if (reason) { 21352392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 2147a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 215a7cc72afSBarry Smith breakout = PETSC_TRUE; 216454a90a3SBarry Smith break; 21752392280SLois Curfman McInnes } 21850ffb88aSMatthew Knepley snes->numFailures++; 2194800dd8cSBarry Smith } 2201acabf8cSLois Curfman McInnes if (!breakout) { 22185385478SLisandro Dalcin /* Update function and solution vectors */ 2224800dd8cSBarry Smith fnorm = gnorm; 22385385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 22485385478SLisandro Dalcin ierr = VecCopy(Y,X);CHKERRQ(ierr); 22585385478SLisandro Dalcin /* Monitor convergence */ 226e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 227c509a162SBarry Smith snes->iter = i+1; 228fbe28522SBarry Smith snes->norm = fnorm; 229e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 230a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2317a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 23285385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 233ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 234e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 235e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 2363f149594SLisandro Dalcin if (reason) break; 2371aa26658SKarl Rupp } else break; 23838442cffSBarry Smith } 23952392280SLois Curfman McInnes if (i == maxits) { 240ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 24185385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 24252392280SLois Curfman McInnes } 243e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 244184914b5SBarry Smith snes->reason = reason; 245e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2474800dd8cSBarry Smith } 2484800dd8cSBarry Smith /*------------------------------------------------------------*/ 2494a2ae208SSatish Balay #undef __FUNCT__ 25004d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONTR" 25104d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 2524800dd8cSBarry Smith { 253dfbe8321SBarry Smith PetscErrorCode ierr; 2543a40ed3dSBarry Smith 2553a40ed3dSBarry Smith PetscFunctionBegin; 256fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 2576cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2583a40ed3dSBarry Smith PetscFunctionReturn(0); 2594800dd8cSBarry Smith } 2606b8b9a38SLisandro Dalcin 2616b8b9a38SLisandro Dalcin #undef __FUNCT__ 26204d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONTR" 26304d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 2646b8b9a38SLisandro Dalcin { 2656b8b9a38SLisandro Dalcin 2666b8b9a38SLisandro Dalcin PetscFunctionBegin; 2676b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2686b8b9a38SLisandro Dalcin } 2696b8b9a38SLisandro Dalcin 2704a2ae208SSatish Balay #undef __FUNCT__ 27104d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONTR" 27204d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 2734800dd8cSBarry Smith { 274dfbe8321SBarry Smith PetscErrorCode ierr; 2755baf8537SBarry Smith 2763a40ed3dSBarry Smith PetscFunctionBegin; 27704d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 278606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 2804800dd8cSBarry Smith } 2814800dd8cSBarry Smith /*------------------------------------------------------------*/ 2824800dd8cSBarry Smith 2834a2ae208SSatish Balay #undef __FUNCT__ 28404d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONTR" 2858c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptions *PetscOptionsObject,SNES snes) 2864800dd8cSBarry Smith { 28704d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 288dfbe8321SBarry Smith PetscErrorCode ierr; 2894800dd8cSBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 291e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 29294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 29394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 29494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 29594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 29694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 29794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 29894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 29994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 300b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3013a40ed3dSBarry Smith PetscFunctionReturn(0); 3024800dd8cSBarry Smith } 3034800dd8cSBarry Smith 3044a2ae208SSatish Balay #undef __FUNCT__ 30504d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONTR" 30604d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 307a935fc98SLois Curfman McInnes { 30804d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 309dfbe8321SBarry Smith PetscErrorCode ierr; 310ace3abfcSBarry Smith PetscBool iascii; 311a935fc98SLois Curfman McInnes 3123a40ed3dSBarry Smith PetscFunctionBegin; 313251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 31432077d6dSBarry Smith if (iascii) { 31557622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 31657622a8eSBarry 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); 31719bcc07fSBarry Smith } 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 319a935fc98SLois Curfman McInnes } 32052392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 321ebe3b25bSBarry Smith /*MC 32204d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 323ebe3b25bSBarry Smith 324ebe3b25bSBarry Smith Options Database: 325ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 326ebe3b25bSBarry Smith . -snes_tr_mu <mu> 327ebe3b25bSBarry Smith . -snes_tr_eta <eta> 328ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 329ebe3b25bSBarry Smith . -snes_tr_delta0 <delta0> 330ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 331ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 332ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 333ebe3b25bSBarry Smith 334ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 335ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 336ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 337ebe3b25bSBarry Smith 338ebe3b25bSBarry Smith This is intended as a model implementation, since it does not 339ebe3b25bSBarry Smith necessarily have many of the bells and whistles of other 340ebe3b25bSBarry Smith implementations. 341ebe3b25bSBarry Smith 342ee3001cbSBarry Smith Level: intermediate 343ee3001cbSBarry Smith 34404d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 345ebe3b25bSBarry Smith 346ebe3b25bSBarry Smith M*/ 3474a2ae208SSatish Balay #undef __FUNCT__ 34804d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONTR" 3498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 3504800dd8cSBarry Smith { 35104d7464bSBarry Smith SNES_NEWTONTR *neP; 352dfbe8321SBarry Smith PetscErrorCode ierr; 3534800dd8cSBarry Smith 3543a40ed3dSBarry Smith PetscFunctionBegin; 35504d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 35604d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 35704d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 35804d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 35904d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 36004d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 361fbe28522SBarry Smith 36242f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 36342f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 36442f4f86dSBarry Smith 365b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 366fbe28522SBarry Smith snes->data = (void*)neP; 367fbe28522SBarry Smith neP->mu = 0.25; 368fbe28522SBarry Smith neP->eta = 0.75; 369fbe28522SBarry Smith neP->delta = 0.0; 370fbe28522SBarry Smith neP->delta0 = 0.2; 371fbe28522SBarry Smith neP->delta1 = 0.3; 372fbe28522SBarry Smith neP->delta2 = 0.75; 373fbe28522SBarry Smith neP->delta3 = 2.0; 374fbe28522SBarry Smith neP->sigma = 0.0001; 375ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 37675567043SBarry Smith neP->rnorm0 = 0.0; 37775567043SBarry Smith neP->ttol = 0.0; 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 3794800dd8cSBarry Smith } 38082bf6240SBarry Smith 381