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 */ 13470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 14df60cc22SBarry Smith { 15971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 16971273eeSBarry Smith SNES snes = ctx->snes; 1704d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 18df60cc22SBarry Smith Vec x; 19064f8208SBarry Smith PetscReal nrm; 20dfbe8321SBarry Smith PetscErrorCode ierr; 21df60cc22SBarry Smith 223a40ed3dSBarry Smith PetscFunctionBegin; 238de4e6dcSJed Brown ierr = KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr); 24329f5518SBarry Smith if (*reason) { 2557622a8eSBarry Smith ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 269b04593eSLois Curfman McInnes } 27a935fc98SLois Curfman McInnes /* Determine norm of solution */ 2878b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 29064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 30064f8208SBarry Smith if (nrm >= neP->delta) { 3157622a8eSBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 32329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 33df60cc22SBarry Smith } 343a40ed3dSBarry Smith PetscFunctionReturn(0); 35df60cc22SBarry Smith } 3682bf6240SBarry Smith 37470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 38971273eeSBarry Smith { 39971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 40971273eeSBarry Smith PetscErrorCode ierr; 41971273eeSBarry Smith 42971273eeSBarry Smith PetscFunctionBegin; 438de4e6dcSJed Brown ierr = KSPConvergedDefaultDestroy(ctx->ctx);CHKERRQ(ierr); 44971273eeSBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 45971273eeSBarry Smith PetscFunctionReturn(0); 46971273eeSBarry Smith } 47971273eeSBarry Smith 4885385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 4985385478SLisandro Dalcin /* 50470880abSPatrick Sanan SNESTR_Converged_Private -test convergence JUST for 5185385478SLisandro Dalcin the trust region tolerance. 5285385478SLisandro Dalcin 5385385478SLisandro Dalcin */ 54470880abSPatrick Sanan static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 5585385478SLisandro Dalcin { 5604d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 5785385478SLisandro Dalcin PetscErrorCode ierr; 5885385478SLisandro Dalcin 5985385478SLisandro Dalcin PetscFunctionBegin; 6085385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 6185385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 6257622a8eSBarry Smith ierr = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 6385385478SLisandro Dalcin *reason = SNES_CONVERGED_TR_DELTA; 64e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 6585385478SLisandro Dalcin ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 6685385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 6785385478SLisandro Dalcin } 6885385478SLisandro Dalcin PetscFunctionReturn(0); 6985385478SLisandro Dalcin } 7085385478SLisandro Dalcin 7185385478SLisandro Dalcin 72df60cc22SBarry Smith /* 7304d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 74ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 754800dd8cSBarry Smith 76ddbbdb52SLois Curfman McInnes 774800dd8cSBarry Smith */ 7804d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 794800dd8cSBarry Smith { 8004d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 8185385478SLisandro Dalcin Vec X,F,Y,G,Ytmp; 826849ba73SBarry Smith PetscErrorCode ierr; 83a7cc72afSBarry Smith PetscInt maxits,i,lits; 8485385478SLisandro Dalcin PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 8579f36460SBarry Smith PetscScalar cnorm; 86df60cc22SBarry Smith KSP ksp; 873f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 88ace3abfcSBarry Smith PetscBool conv = PETSC_FALSE,breakout = PETSC_FALSE; 894800dd8cSBarry Smith 903a40ed3dSBarry Smith PetscFunctionBegin; 916c4ed002SBarry 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); 92c579b300SPatrick Farrell 93fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 94fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 9539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 96fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 97fbe28522SBarry Smith G = snes->work[1]; 986b5873e3SBarry Smith Ytmp = snes->work[2]; 994800dd8cSBarry Smith 100e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1014c49b128SBarry Smith snes->iter = 0; 102e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1034800dd8cSBarry Smith 104e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1055334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 1061aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 107e4ed7901SPeter Brune 108cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 109422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 1108b84755bSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 111e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 112fbe28522SBarry Smith snes->norm = fnorm; 113e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1148b84755bSBarry Smith delta = xnorm ? neP->delta0*xnorm : neP->delta0; 1154800dd8cSBarry Smith neP->delta = delta; 116a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1177a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 118b37302c6SLois Curfman McInnes 11985385478SLisandro Dalcin /* test convergence */ 12085385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 12185385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 1223f149594SLisandro Dalcin 123a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 124c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 1257c4af3dcSLois Curfman McInnes if (!conv) { 126971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 1273f149594SLisandro Dalcin ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 128b00a9115SJed Brown ierr = PetscNew(&ctx);CHKERRQ(ierr); 129971273eeSBarry Smith ctx->snes = snes; 1308de4e6dcSJed Brown ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 131470880abSPatrick Sanan ierr = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr); 132470880abSPatrick Sanan ierr = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr); 1337c4af3dcSLois Curfman McInnes } 134df60cc22SBarry Smith 1354800dd8cSBarry Smith for (i=0; i<maxits; i++) { 13699a96b7cSMatthew Knepley 13799a96b7cSMatthew Knepley /* Call general purpose update function */ 138e7788613SBarry Smith if (snes->ops->update) { 139e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 14099a96b7cSMatthew Knepley } 14199a96b7cSMatthew Knepley 14285385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 143d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 14407b62357SFande Kong SNESCheckJacobianDomainerror(snes); 14523ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 146d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 1473f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1481aa26658SKarl Rupp 149329f5518SBarry Smith snes->linear_its += lits; 1501aa26658SKarl Rupp 151ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 152064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 153064f8208SBarry Smith norm1 = nrm; 1546b5873e3SBarry Smith while (1) { 155393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 156064f8208SBarry Smith nrm = norm1; 1576b5873e3SBarry Smith 1582deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 159064f8208SBarry Smith if (nrm >= delta) { 160064f8208SBarry Smith nrm = delta/nrm; 161064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 162064f8208SBarry Smith cnorm = nrm; 16357622a8eSBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 1642dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 165064f8208SBarry Smith nrm = gpnorm; 166fbe28522SBarry Smith ynorm = delta; 167fbe28522SBarry Smith } else { 168fbe28522SBarry Smith gpnorm = 0.0; 169ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 170064f8208SBarry Smith ynorm = nrm; 171fbe28522SBarry Smith } 172*1592aa04SStefano Zampini ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 17379f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 1745334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 175cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1764800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 177fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1784800dd8cSBarry Smith 1794800dd8cSBarry Smith /* Update size of trust region */ 1804800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1814800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1824800dd8cSBarry Smith else delta *= neP->delta3; 18357622a8eSBarry Smith ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 18457622a8eSBarry Smith ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 1851aa26658SKarl Rupp 1864800dd8cSBarry Smith neP->delta = delta; 1874800dd8cSBarry Smith if (rho > neP->sigma) break; 188ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 1896b5873e3SBarry Smith /* check to see if progress is hopeless */ 190ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 191470880abSPatrick Sanan ierr = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 19285385478SLisandro Dalcin if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 193184914b5SBarry Smith if (reason) { 19452392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 1957a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 196a7cc72afSBarry Smith breakout = PETSC_TRUE; 197454a90a3SBarry Smith break; 19852392280SLois Curfman McInnes } 19950ffb88aSMatthew Knepley snes->numFailures++; 2004800dd8cSBarry Smith } 2011acabf8cSLois Curfman McInnes if (!breakout) { 20285385478SLisandro Dalcin /* Update function and solution vectors */ 2034800dd8cSBarry Smith fnorm = gnorm; 20485385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 20585385478SLisandro Dalcin ierr = VecCopy(Y,X);CHKERRQ(ierr); 20685385478SLisandro Dalcin /* Monitor convergence */ 207e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 208c509a162SBarry Smith snes->iter = i+1; 209fbe28522SBarry Smith snes->norm = fnorm; 210c1e67a49SFande Kong snes->xnorm = xnorm; 211c1e67a49SFande Kong snes->ynorm = ynorm; 212e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 213a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2147a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 21585385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 216ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 217e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 218e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 2193f149594SLisandro Dalcin if (reason) break; 2201aa26658SKarl Rupp } else break; 22138442cffSBarry Smith } 22252392280SLois Curfman McInnes if (i == maxits) { 223ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 22485385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 22552392280SLois Curfman McInnes } 226e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 227184914b5SBarry Smith snes->reason = reason; 228e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 2304800dd8cSBarry Smith } 2314800dd8cSBarry Smith /*------------------------------------------------------------*/ 23204d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 2334800dd8cSBarry Smith { 234dfbe8321SBarry Smith PetscErrorCode ierr; 2353a40ed3dSBarry Smith 2363a40ed3dSBarry Smith PetscFunctionBegin; 237fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 2386cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2393a40ed3dSBarry Smith PetscFunctionReturn(0); 2404800dd8cSBarry Smith } 2416b8b9a38SLisandro Dalcin 24204d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes) 2436b8b9a38SLisandro Dalcin { 2446b8b9a38SLisandro Dalcin 2456b8b9a38SLisandro Dalcin PetscFunctionBegin; 2466b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2476b8b9a38SLisandro Dalcin } 2486b8b9a38SLisandro Dalcin 24904d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 2504800dd8cSBarry Smith { 251dfbe8321SBarry Smith PetscErrorCode ierr; 2525baf8537SBarry Smith 2533a40ed3dSBarry Smith PetscFunctionBegin; 25404d7464bSBarry Smith ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 255606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2563a40ed3dSBarry Smith PetscFunctionReturn(0); 2574800dd8cSBarry Smith } 2584800dd8cSBarry Smith /*------------------------------------------------------------*/ 2594800dd8cSBarry Smith 2604416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes) 2614800dd8cSBarry Smith { 26204d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 263dfbe8321SBarry Smith PetscErrorCode ierr; 2644800dd8cSBarry Smith 2653a40ed3dSBarry Smith PetscFunctionBegin; 266e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 26794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 26894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 26994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 27094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 27194ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 27294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 27394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 27494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 275b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2763a40ed3dSBarry Smith PetscFunctionReturn(0); 2774800dd8cSBarry Smith } 2784800dd8cSBarry Smith 27904d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 280a935fc98SLois Curfman McInnes { 28104d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 282dfbe8321SBarry Smith PetscErrorCode ierr; 283ace3abfcSBarry Smith PetscBool iascii; 284a935fc98SLois Curfman McInnes 2853a40ed3dSBarry Smith PetscFunctionBegin; 286251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 28732077d6dSBarry Smith if (iascii) { 28854fe5c21SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr); 28957622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 29057622a8eSBarry 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); 29119bcc07fSBarry Smith } 2923a40ed3dSBarry Smith PetscFunctionReturn(0); 293a935fc98SLois Curfman McInnes } 29452392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 295ebe3b25bSBarry Smith /*MC 29604d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 297ebe3b25bSBarry Smith 298ebe3b25bSBarry Smith Options Database: 2998e434772SBarry Smith + -snes_trtol <tol> - trust region tolerance 3008e434772SBarry Smith . -snes_tr_mu <mu> - trust region parameter 3018e434772SBarry Smith . -snes_tr_eta <eta> - trust region parameter 3028e434772SBarry Smith . -snes_tr_sigma <sigma> - trust region parameter 3038b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 3048e434772SBarry Smith . -snes_tr_delta1 <delta1> - trust region parameter 3058e434772SBarry Smith . -snes_tr_delta2 <delta2> - trust region parameter 3068e434772SBarry Smith - -snes_tr_delta3 <delta3> - trust region parameter 307ebe3b25bSBarry Smith 308ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 309ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 310ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 311ebe3b25bSBarry Smith 312ee3001cbSBarry Smith Level: intermediate 313ee3001cbSBarry Smith 31404d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 315ebe3b25bSBarry Smith 316ebe3b25bSBarry Smith M*/ 3178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 3184800dd8cSBarry Smith { 31904d7464bSBarry Smith SNES_NEWTONTR *neP; 320dfbe8321SBarry Smith PetscErrorCode ierr; 3214800dd8cSBarry Smith 3223a40ed3dSBarry Smith PetscFunctionBegin; 32304d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 32404d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 32504d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 32604d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 32704d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 32804d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 329fbe28522SBarry Smith 33042f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 331efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 33242f4f86dSBarry Smith 3334fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3344fc747eaSLawrence Mitchell 335b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 336fbe28522SBarry Smith snes->data = (void*)neP; 337fbe28522SBarry Smith neP->mu = 0.25; 338fbe28522SBarry Smith neP->eta = 0.75; 339fbe28522SBarry Smith neP->delta = 0.0; 340fbe28522SBarry Smith neP->delta0 = 0.2; 341fbe28522SBarry Smith neP->delta1 = 0.3; 342fbe28522SBarry Smith neP->delta2 = 0.75; 343fbe28522SBarry Smith neP->delta3 = 2.0; 344fbe28522SBarry Smith neP->sigma = 0.0001; 345ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 34675567043SBarry Smith neP->rnorm0 = 0.0; 34775567043SBarry Smith neP->ttol = 0.0; 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 3494800dd8cSBarry Smith } 35082bf6240SBarry Smith 351