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; 25*8de4e6dcSJed Brown ierr = KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr); 26329f5518SBarry Smith if (*reason) { 2771f87433Sdalcinl ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,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) { 33ae15b995SBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,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; 47*8de4e6dcSJed 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) { 6885385478SLisandro Dalcin ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,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; 92112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 9385385478SLisandro Dalcin PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 9479f36460SBarry Smith PetscScalar cnorm; 95df60cc22SBarry Smith KSP ksp; 963f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 97ace3abfcSBarry Smith PetscBool conv = PETSC_FALSE,breakout = PETSC_FALSE; 98e4ed7901SPeter Brune PetscBool domainerror; 994800dd8cSBarry Smith 1003a40ed3dSBarry Smith PetscFunctionBegin; 101fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 102fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 10339e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 104fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 105fbe28522SBarry Smith G = snes->work[1]; 1066b5873e3SBarry Smith Ytmp = snes->work[2]; 1074800dd8cSBarry Smith 108ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1094c49b128SBarry Smith snes->iter = 0; 110ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1114800dd8cSBarry Smith 112e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1135334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 114e4ed7901SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 115e4ed7901SPeter Brune if (domainerror) { 116e4ed7901SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 117e4ed7901SPeter Brune PetscFunctionReturn(0); 118e4ed7901SPeter Brune } 1191aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 120e4ed7901SPeter Brune 121e4ed7901SPeter Brune if (!snes->norm_init_set) { 122cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 123189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 124189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 125189a9710SBarry Smith PetscFunctionReturn(0); 126189a9710SBarry Smith } 127e4ed7901SPeter Brune } else { 128e4ed7901SPeter Brune fnorm = snes->norm_init; 129e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 130e4ed7901SPeter Brune } 131e4ed7901SPeter Brune 132ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 133fbe28522SBarry Smith snes->norm = fnorm; 134ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1354800dd8cSBarry Smith delta = neP->delta0*fnorm; 1364800dd8cSBarry Smith neP->delta = delta; 137a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1387a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 139b37302c6SLois Curfman McInnes 14085385478SLisandro Dalcin /* test convergence */ 14185385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 14285385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 1433f149594SLisandro Dalcin 144a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 1450298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 1467c4af3dcSLois Curfman McInnes if (!conv) { 147971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 1483f149594SLisandro Dalcin ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 149971273eeSBarry Smith ierr = PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);CHKERRQ(ierr); 150971273eeSBarry Smith ctx->snes = snes; 151*8de4e6dcSJed Brown ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 152971273eeSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr); 153ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 1547c4af3dcSLois Curfman McInnes } 155df60cc22SBarry Smith 1564800dd8cSBarry Smith for (i=0; i<maxits; i++) { 15799a96b7cSMatthew Knepley 15899a96b7cSMatthew Knepley /* Call general purpose update function */ 159e7788613SBarry Smith if (snes->ops->update) { 160e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 16199a96b7cSMatthew Knepley } 16299a96b7cSMatthew Knepley 16385385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 1645334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 16594b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 166d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 1673f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1681aa26658SKarl Rupp 169329f5518SBarry Smith snes->linear_its += lits; 1701aa26658SKarl Rupp 171ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 172064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 173064f8208SBarry Smith norm1 = nrm; 1746b5873e3SBarry Smith while (1) { 175393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 176064f8208SBarry Smith nrm = norm1; 1776b5873e3SBarry Smith 1782deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 179064f8208SBarry Smith if (nrm >= delta) { 180064f8208SBarry Smith nrm = delta/nrm; 181064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 182064f8208SBarry Smith cnorm = nrm; 183ae15b995SBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr); 1842dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 185064f8208SBarry Smith nrm = gpnorm; 186fbe28522SBarry Smith ynorm = delta; 187fbe28522SBarry Smith } else { 188fbe28522SBarry Smith gpnorm = 0.0; 189ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 190064f8208SBarry Smith ynorm = nrm; 191fbe28522SBarry Smith } 19279f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 19385385478SLisandro Dalcin ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 1945334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 195cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1964800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 197fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1984800dd8cSBarry Smith 1994800dd8cSBarry Smith /* Update size of trust region */ 2004800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 2014800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 2024800dd8cSBarry Smith else delta *= neP->delta3; 203ae15b995SBarry Smith ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr); 204ae15b995SBarry Smith ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr); 2051aa26658SKarl Rupp 2064800dd8cSBarry Smith neP->delta = delta; 2074800dd8cSBarry Smith if (rho > neP->sigma) break; 208ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 2096b5873e3SBarry Smith /* check to see if progress is hopeless */ 210ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 21185385478SLisandro Dalcin ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 21285385478SLisandro Dalcin if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 213184914b5SBarry Smith if (reason) { 21452392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 2157a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 216a7cc72afSBarry Smith breakout = PETSC_TRUE; 217454a90a3SBarry Smith break; 21852392280SLois Curfman McInnes } 21950ffb88aSMatthew Knepley snes->numFailures++; 2204800dd8cSBarry Smith } 2211acabf8cSLois Curfman McInnes if (!breakout) { 22285385478SLisandro Dalcin /* Update function and solution vectors */ 2234800dd8cSBarry Smith fnorm = gnorm; 22485385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 22585385478SLisandro Dalcin ierr = VecCopy(Y,X);CHKERRQ(ierr); 22685385478SLisandro Dalcin /* Monitor convergence */ 227ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 228c509a162SBarry Smith snes->iter = i+1; 229fbe28522SBarry Smith snes->norm = fnorm; 230ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 231a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2327a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 23385385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 234ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 235e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 236e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 2373f149594SLisandro Dalcin if (reason) break; 2381aa26658SKarl Rupp } else break; 23938442cffSBarry Smith } 24052392280SLois Curfman McInnes if (i == maxits) { 241ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 24285385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 24352392280SLois Curfman McInnes } 244ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 245184914b5SBarry Smith snes->reason = reason; 246ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2473a40ed3dSBarry Smith PetscFunctionReturn(0); 2484800dd8cSBarry Smith } 2494800dd8cSBarry Smith /*------------------------------------------------------------*/ 2504a2ae208SSatish Balay #undef __FUNCT__ 25104d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONTR" 25204d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 2534800dd8cSBarry Smith { 254dfbe8321SBarry Smith PetscErrorCode ierr; 2553a40ed3dSBarry Smith 2563a40ed3dSBarry Smith PetscFunctionBegin; 257fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 2586cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2593a40ed3dSBarry Smith PetscFunctionReturn(0); 2604800dd8cSBarry Smith } 2616b8b9a38SLisandro Dalcin 2626b8b9a38SLisandro Dalcin #undef __FUNCT__ 26304d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONTR" 26404d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 2656b8b9a38SLisandro Dalcin { 2666b8b9a38SLisandro Dalcin 2676b8b9a38SLisandro Dalcin PetscFunctionBegin; 2686b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2696b8b9a38SLisandro Dalcin } 2706b8b9a38SLisandro Dalcin 2714a2ae208SSatish Balay #undef __FUNCT__ 27204d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONTR" 27304d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 2744800dd8cSBarry Smith { 275dfbe8321SBarry Smith PetscErrorCode ierr; 2765baf8537SBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 27804d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 279606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2803a40ed3dSBarry Smith PetscFunctionReturn(0); 2814800dd8cSBarry Smith } 2824800dd8cSBarry Smith /*------------------------------------------------------------*/ 2834800dd8cSBarry Smith 2844a2ae208SSatish Balay #undef __FUNCT__ 28504d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONTR" 28604d7464bSBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes) 2874800dd8cSBarry Smith { 28804d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 289dfbe8321SBarry Smith PetscErrorCode ierr; 2904800dd8cSBarry Smith 2913a40ed3dSBarry Smith PetscFunctionBegin; 292b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 29387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 2944b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 2954b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 2964b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 2974b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 2984b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 2994b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 3004b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 301b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 3034800dd8cSBarry Smith } 3044800dd8cSBarry Smith 3054a2ae208SSatish Balay #undef __FUNCT__ 30604d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONTR" 30704d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 308a935fc98SLois Curfman McInnes { 30904d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 310dfbe8321SBarry Smith PetscErrorCode ierr; 311ace3abfcSBarry Smith PetscBool iascii; 312a935fc98SLois Curfman McInnes 3133a40ed3dSBarry Smith PetscFunctionBegin; 314251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 31532077d6dSBarry Smith if (iascii) { 316a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 317a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 31819bcc07fSBarry Smith } 3193a40ed3dSBarry Smith PetscFunctionReturn(0); 320a935fc98SLois Curfman McInnes } 32152392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 322ebe3b25bSBarry Smith /*MC 32304d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 324ebe3b25bSBarry Smith 325ebe3b25bSBarry Smith Options Database: 326ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 327ebe3b25bSBarry Smith . -snes_tr_mu <mu> 328ebe3b25bSBarry Smith . -snes_tr_eta <eta> 329ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 330ebe3b25bSBarry Smith . -snes_tr_delta0 <delta0> 331ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 332ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 333ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 334ebe3b25bSBarry Smith 335ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 336ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 337ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 338ebe3b25bSBarry Smith 339ebe3b25bSBarry Smith This is intended as a model implementation, since it does not 340ebe3b25bSBarry Smith necessarily have many of the bells and whistles of other 341ebe3b25bSBarry Smith implementations. 342ebe3b25bSBarry Smith 343ee3001cbSBarry Smith Level: intermediate 344ee3001cbSBarry Smith 34504d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 346ebe3b25bSBarry Smith 347ebe3b25bSBarry Smith M*/ 3484a2ae208SSatish Balay #undef __FUNCT__ 34904d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONTR" 3508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 3514800dd8cSBarry Smith { 35204d7464bSBarry Smith SNES_NEWTONTR *neP; 353dfbe8321SBarry Smith PetscErrorCode ierr; 3544800dd8cSBarry Smith 3553a40ed3dSBarry Smith PetscFunctionBegin; 35604d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 35704d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 35804d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 35904d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 36004d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 36104d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 362fbe28522SBarry Smith 36342f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 36442f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 36542f4f86dSBarry Smith 36604d7464bSBarry Smith ierr = PetscNewLog(snes,SNES_NEWTONTR,&neP);CHKERRQ(ierr); 367fbe28522SBarry Smith snes->data = (void*)neP; 368fbe28522SBarry Smith neP->mu = 0.25; 369fbe28522SBarry Smith neP->eta = 0.75; 370fbe28522SBarry Smith neP->delta = 0.0; 371fbe28522SBarry Smith neP->delta0 = 0.2; 372fbe28522SBarry Smith neP->delta1 = 0.3; 373fbe28522SBarry Smith neP->delta2 = 0.75; 374fbe28522SBarry Smith neP->delta3 = 2.0; 375fbe28522SBarry Smith neP->sigma = 0.0001; 376ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 37775567043SBarry Smith neP->rnorm0 = 0.0; 37875567043SBarry Smith neP->ttol = 0.0; 3793a40ed3dSBarry Smith PetscFunctionReturn(0); 3804800dd8cSBarry Smith } 38182bf6240SBarry Smith 382