163dd3a1aSKris Buschelman #define PETSCSNES_DLL 24800dd8cSBarry Smith 3e090d566SSatish Balay #include "src/snes/impls/tr/tr.h" /*I "petscsnes.h" I*/ 44800dd8cSBarry Smith 5fbe28522SBarry Smith /* 6df60cc22SBarry Smith This convergence test determines if the two norm of the 7df60cc22SBarry Smith solution lies outside the trust region, if so it halts. 8df60cc22SBarry Smith */ 94a2ae208SSatish Balay #undef __FUNCT__ 104b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private" 1177431f27SBarry Smith PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx) 12df60cc22SBarry Smith { 13df60cc22SBarry Smith SNES snes = (SNES) ctx; 144b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR*)snes->data; 15df60cc22SBarry Smith Vec x; 16064f8208SBarry Smith PetscReal nrm; 17dfbe8321SBarry Smith PetscErrorCode ierr; 18df60cc22SBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 20329f5518SBarry Smith ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr); 21329f5518SBarry Smith if (*reason) { 2271f87433Sdalcinl ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr); 239b04593eSLois Curfman McInnes } 24a935fc98SLois Curfman McInnes /* Determine norm of solution */ 2578b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 26064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 27064f8208SBarry Smith if (nrm >= neP->delta) { 28ae15b995SBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr); 29329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 30df60cc22SBarry Smith } 313a40ed3dSBarry Smith PetscFunctionReturn(0); 32df60cc22SBarry Smith } 3382bf6240SBarry Smith 34*85385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 35*85385478SLisandro Dalcin #undef __FUNCT__ 36*85385478SLisandro Dalcin #define __FUNCT__ "SNES_TR_Converged_Private" 37*85385478SLisandro Dalcin /* 38*85385478SLisandro Dalcin SNES_TR_Converged_Private -test convergence JUST for 39*85385478SLisandro Dalcin the trust region tolerance. 40*85385478SLisandro Dalcin 41*85385478SLisandro Dalcin */ 42*85385478SLisandro Dalcin static PetscErrorCode SNES_TR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 43*85385478SLisandro Dalcin { 44*85385478SLisandro Dalcin SNES_TR *neP = (SNES_TR *)snes->data; 45*85385478SLisandro Dalcin PetscErrorCode ierr; 46*85385478SLisandro Dalcin 47*85385478SLisandro Dalcin PetscFunctionBegin; 48*85385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 49*85385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 50*85385478SLisandro Dalcin ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr); 51*85385478SLisandro Dalcin *reason = SNES_CONVERGED_TR_DELTA; 52*85385478SLisandro Dalcin } else if (snes->nfuncs >= snes->max_funcs) { 53*85385478SLisandro Dalcin ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 54*85385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 55*85385478SLisandro Dalcin } 56*85385478SLisandro Dalcin PetscFunctionReturn(0); 57*85385478SLisandro Dalcin } 58*85385478SLisandro Dalcin 59*85385478SLisandro Dalcin 60df60cc22SBarry Smith /* 614b27c08aSLois Curfman McInnes SNESSolve_TR - Implements Newton's Method with a very simple trust 62ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 634800dd8cSBarry Smith 64ddbbdb52SLois Curfman McInnes 654800dd8cSBarry Smith */ 664a2ae208SSatish Balay #undef __FUNCT__ 674b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR" 686849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes) 694800dd8cSBarry Smith { 704b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR*)snes->data; 71*85385478SLisandro Dalcin Vec X,F,Y,G,Ytmp; 726849ba73SBarry Smith PetscErrorCode ierr; 73a7cc72afSBarry Smith PetscInt maxits,i,lits; 74112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 75*85385478SLisandro Dalcin PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 7679f36460SBarry Smith PetscScalar cnorm; 77df60cc22SBarry Smith KSP ksp; 783f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 79a7cc72afSBarry Smith PetscTruth conv,breakout = PETSC_FALSE; 804800dd8cSBarry Smith 813a40ed3dSBarry Smith PetscFunctionBegin; 82fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 83fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 8439e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 85fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 86fbe28522SBarry Smith G = snes->work[1]; 876b5873e3SBarry Smith Ytmp = snes->work[2]; 884800dd8cSBarry Smith 894c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 904c49b128SBarry Smith snes->iter = 0; 914c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 924800dd8cSBarry Smith 935334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 94cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 950f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 96fbe28522SBarry Smith snes->norm = fnorm; 970f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 984800dd8cSBarry Smith delta = neP->delta0*fnorm; 994800dd8cSBarry Smith neP->delta = delta; 1000462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 10194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 102b37302c6SLois Curfman McInnes 103d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 104d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 105*85385478SLisandro Dalcin /* test convergence */ 106*85385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 107*85385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 1083f149594SLisandro Dalcin 109a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 1107c4af3dcSLois Curfman McInnes ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 1117c4af3dcSLois Curfman McInnes if (!conv) { 1123f149594SLisandro Dalcin ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1134b27c08aSLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr); 114ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 1157c4af3dcSLois Curfman McInnes } 116df60cc22SBarry Smith 1174800dd8cSBarry Smith for (i=0; i<maxits; i++) { 11899a96b7cSMatthew Knepley 11999a96b7cSMatthew Knepley /* Call general purpose update function */ 120e7788613SBarry Smith if (snes->ops->update) { 121e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 12299a96b7cSMatthew Knepley } 12399a96b7cSMatthew Knepley 124*85385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 1255334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 12694b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1273f149594SLisandro Dalcin ierr = SNES_KSPSolve(snes,snes->ksp,F,Ytmp);CHKERRQ(ierr); 1283f149594SLisandro Dalcin ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 129329f5518SBarry Smith snes->linear_its += lits; 130ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 131064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 132064f8208SBarry Smith norm1 = nrm; 1336b5873e3SBarry Smith while(1) { 134393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 135064f8208SBarry Smith nrm = norm1; 1366b5873e3SBarry Smith 1372deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 138064f8208SBarry Smith if (nrm >= delta) { 139064f8208SBarry Smith nrm = delta/nrm; 140064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 141064f8208SBarry Smith cnorm = nrm; 142ae15b995SBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr); 1432dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 144064f8208SBarry Smith nrm = gpnorm; 145fbe28522SBarry Smith ynorm = delta; 146fbe28522SBarry Smith } else { 147fbe28522SBarry Smith gpnorm = 0.0; 148ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 149064f8208SBarry Smith ynorm = nrm; 150fbe28522SBarry Smith } 15179f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 152*85385478SLisandro Dalcin ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 1535334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 154cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1554800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 156fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1574800dd8cSBarry Smith 1584800dd8cSBarry Smith /* Update size of trust region */ 1594800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1604800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1614800dd8cSBarry Smith else delta *= neP->delta3; 162ae15b995SBarry Smith ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr); 163ae15b995SBarry Smith ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr); 1644800dd8cSBarry Smith neP->delta = delta; 1654800dd8cSBarry Smith if (rho > neP->sigma) break; 166ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 1676b5873e3SBarry Smith /* check to see if progress is hopeless */ 168ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 169*85385478SLisandro Dalcin ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 170*85385478SLisandro Dalcin if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 171184914b5SBarry Smith if (reason) { 17252392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 1735ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 174a7cc72afSBarry Smith breakout = PETSC_TRUE; 175454a90a3SBarry Smith break; 17652392280SLois Curfman McInnes } 17750ffb88aSMatthew Knepley snes->numFailures++; 1784800dd8cSBarry Smith } 1791acabf8cSLois Curfman McInnes if (!breakout) { 180*85385478SLisandro Dalcin /* Update function and solution vectors */ 1814800dd8cSBarry Smith fnorm = gnorm; 182*85385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 183*85385478SLisandro Dalcin ierr = VecCopy(Y,X);CHKERRQ(ierr); 184*85385478SLisandro Dalcin /* Monitor convergence */ 1850f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 186c509a162SBarry Smith snes->iter = i+1; 187fbe28522SBarry Smith snes->norm = fnorm; 1880f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 189*85385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 190*85385478SLisandro Dalcin SNESMonitor(snes,snes->iter,snes->norm); 191*85385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 192ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 193*85385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 194e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 1953f149594SLisandro Dalcin if (reason) break; 1961acabf8cSLois Curfman McInnes } else { 1971acabf8cSLois Curfman McInnes break; 1981acabf8cSLois Curfman McInnes } 19938442cffSBarry Smith } 20052392280SLois Curfman McInnes if (i == maxits) { 201ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 202*85385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 20352392280SLois Curfman McInnes } 2040f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 205184914b5SBarry Smith snes->reason = reason; 2060f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2073a40ed3dSBarry Smith PetscFunctionReturn(0); 2084800dd8cSBarry Smith } 2094800dd8cSBarry Smith /*------------------------------------------------------------*/ 2104a2ae208SSatish Balay #undef __FUNCT__ 2114b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR" 2126849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes) 2134800dd8cSBarry Smith { 214dfbe8321SBarry Smith PetscErrorCode ierr; 2153a40ed3dSBarry Smith 2163a40ed3dSBarry Smith PetscFunctionBegin; 217*85385478SLisandro Dalcin if (!snes->vec_sol_update) { 218*85385478SLisandro Dalcin ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 219*85385478SLisandro Dalcin ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 220*85385478SLisandro Dalcin } 221410397dcSLisandro Dalcin if (!snes->work) { 222*85385478SLisandro Dalcin snes->nwork = 3; 223d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 224efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 225410397dcSLisandro Dalcin } 2263a40ed3dSBarry Smith PetscFunctionReturn(0); 2274800dd8cSBarry Smith } 2284800dd8cSBarry Smith /*------------------------------------------------------------*/ 2294a2ae208SSatish Balay #undef __FUNCT__ 2304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR" 2316849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes) 2324800dd8cSBarry Smith { 233dfbe8321SBarry Smith PetscErrorCode ierr; 2345baf8537SBarry Smith 2353a40ed3dSBarry Smith PetscFunctionBegin; 236*85385478SLisandro Dalcin if (snes->vec_sol_update) { 237*85385478SLisandro Dalcin ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 238*85385478SLisandro Dalcin snes->vec_sol_update = PETSC_NULL; 239*85385478SLisandro Dalcin } 2405baf8537SBarry Smith if (snes->nwork) { 2414b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 242*85385478SLisandro Dalcin snes->nwork = 0; 243*85385478SLisandro Dalcin snes->work = PETSC_NULL; 2445baf8537SBarry Smith } 245606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2474800dd8cSBarry Smith } 2484800dd8cSBarry Smith /*------------------------------------------------------------*/ 2494800dd8cSBarry Smith 2504a2ae208SSatish Balay #undef __FUNCT__ 2514b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR" 2526849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes) 2534800dd8cSBarry Smith { 2544b27c08aSLois Curfman McInnes SNES_TR *ctx = (SNES_TR *)snes->data; 255dfbe8321SBarry Smith PetscErrorCode ierr; 2564800dd8cSBarry Smith 2573a40ed3dSBarry Smith PetscFunctionBegin; 258b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 25987828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 2604b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 2614b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 2624b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 2634b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 2644b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 2654b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 2664b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 267b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2683a40ed3dSBarry Smith PetscFunctionReturn(0); 2694800dd8cSBarry Smith } 2704800dd8cSBarry Smith 2714a2ae208SSatish Balay #undef __FUNCT__ 2724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR" 2736849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer) 274a935fc98SLois Curfman McInnes { 2754b27c08aSLois Curfman McInnes SNES_TR *tr = (SNES_TR *)snes->data; 276dfbe8321SBarry Smith PetscErrorCode ierr; 27732077d6dSBarry Smith PetscTruth iascii; 278a935fc98SLois Curfman McInnes 2793a40ed3dSBarry Smith PetscFunctionBegin; 28032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 28132077d6dSBarry Smith if (iascii) { 282a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 283a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2845cd90555SBarry Smith } else { 28579a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 28619bcc07fSBarry Smith } 2873a40ed3dSBarry Smith PetscFunctionReturn(0); 288a935fc98SLois Curfman McInnes } 28952392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 290ebe3b25bSBarry Smith /*MC 291ebe3b25bSBarry Smith SNESTR - Newton based nonlinear solver that uses a trust region 292ebe3b25bSBarry Smith 293ebe3b25bSBarry Smith Options Database: 294ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 295ebe3b25bSBarry Smith . -snes_tr_mu <mu> 296ebe3b25bSBarry Smith . -snes_tr_eta <eta> 297ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 298ebe3b25bSBarry Smith . -snes_tr_delta0 <delta0> 299ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 300ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 301ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 302ebe3b25bSBarry Smith 303ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 304ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 305ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 306ebe3b25bSBarry Smith 307ebe3b25bSBarry Smith This is intended as a model implementation, since it does not 308ebe3b25bSBarry Smith necessarily have many of the bells and whistles of other 309ebe3b25bSBarry Smith implementations. 310ebe3b25bSBarry Smith 311ee3001cbSBarry Smith Level: intermediate 312ee3001cbSBarry Smith 313ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance() 314ebe3b25bSBarry Smith 315ebe3b25bSBarry Smith M*/ 316fb2e594dSBarry Smith EXTERN_C_BEGIN 3174a2ae208SSatish Balay #undef __FUNCT__ 3184b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR" 31963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_TR(SNES snes) 3204800dd8cSBarry Smith { 3214b27c08aSLois Curfman McInnes SNES_TR *neP; 322dfbe8321SBarry Smith PetscErrorCode ierr; 3234800dd8cSBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 325e7788613SBarry Smith snes->ops->setup = SNESSetUp_TR; 326e7788613SBarry Smith snes->ops->solve = SNESSolve_TR; 327e7788613SBarry Smith snes->ops->destroy = SNESDestroy_TR; 328e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_TR; 329e7788613SBarry Smith snes->ops->view = SNESView_TR; 330fbe28522SBarry Smith 33138f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_TR,&neP);CHKERRQ(ierr); 332fbe28522SBarry Smith snes->data = (void*)neP; 333fbe28522SBarry Smith neP->mu = 0.25; 334fbe28522SBarry Smith neP->eta = 0.75; 335fbe28522SBarry Smith neP->delta = 0.0; 336fbe28522SBarry Smith neP->delta0 = 0.2; 337fbe28522SBarry Smith neP->delta1 = 0.3; 338fbe28522SBarry Smith neP->delta2 = 0.75; 339fbe28522SBarry Smith neP->delta3 = 2.0; 340fbe28522SBarry Smith neP->sigma = 0.0001; 341ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 34252392280SLois Curfman McInnes neP->rnorm0 = 0; 34352392280SLois Curfman McInnes neP->ttol = 0; 3443a40ed3dSBarry Smith PetscFunctionReturn(0); 3454800dd8cSBarry Smith } 346fb2e594dSBarry Smith EXTERN_C_END 34782bf6240SBarry Smith 348