1*29bbc08cSBarry Smith /*$Id: tr.c,v 1.119 2000/09/02 02:49:40 bsmith Exp bsmith $*/ 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 */ 95615d1e5SSatish Balay #undef __FUNC__ 10b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNES_EQ_TR_KSPConverged_Private" 113acb151aSSatish Balay int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,double rnorm,KSPConvergedReason *reason,void *ctx) 12df60cc22SBarry Smith { 13df60cc22SBarry Smith SNES snes = (SNES) ctx; 1498120f82SLois Curfman McInnes SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 156831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data; 16df60cc22SBarry Smith Vec x; 1798120f82SLois Curfman McInnes double norm; 183acb151aSSatish Balay int ierr; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 2198120f82SLois Curfman McInnes if (snes->ksp_ewconv) { 22*29bbc08cSBarry Smith if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence context not created"); 23c50d6201SSatish Balay if (!n) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);} 2498120f82SLois Curfman McInnes kctx->lresid_last = rnorm; 25df60cc22SBarry Smith } 26329f5518SBarry Smith ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr); 27329f5518SBarry Smith if (*reason) { 28329f5518SBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%g\n",n,rnorm); 299b04593eSLois Curfman McInnes } 30df60cc22SBarry Smith 31a935fc98SLois Curfman McInnes /* Determine norm of solution */ 3278b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 33cddf8d76SBarry Smith ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 34df60cc22SBarry Smith if (norm >= neP->delta) { 356831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 36c38d4ed2SBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,norm); 37329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 38df60cc22SBarry Smith } 393a40ed3dSBarry Smith PetscFunctionReturn(0); 40df60cc22SBarry Smith } 4182bf6240SBarry Smith 42df60cc22SBarry Smith /* 43f63b844aSLois Curfman McInnes SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust 44ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 454800dd8cSBarry Smith 464800dd8cSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 474800dd8cSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 48ddbbdb52SLois Curfman McInnes of Mathematical Software", Wayne Cowell, editor. 49ddbbdb52SLois Curfman McInnes 504800dd8cSBarry Smith This is intended as a model implementation, since it does not 514800dd8cSBarry Smith necessarily have many of the bells and whistles of other 524800dd8cSBarry Smith implementations. 534800dd8cSBarry Smith */ 545615d1e5SSatish Balay #undef __FUNC__ 55b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_TR" 56f63b844aSLois Curfman McInnes static int SNESSolve_EQ_TR(SNES snes,int *its) 574800dd8cSBarry Smith { 586831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data; 596b5873e3SBarry Smith Vec X,F,Y,G,TMP,Ytmp; 600462333dSBarry Smith int maxits,i,ierr,lits,breakout = 0; 61112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 620462333dSBarry Smith double rho,fnorm,gnorm,gpnorm,xnorm,delta,norm,ynorm,norm1; 635334005bSBarry Smith Scalar mone = -1.0,cnorm; 64df60cc22SBarry Smith KSP ksp; 65df60cc22SBarry Smith SLES sles; 66184914b5SBarry Smith SNESConvergedReason reason; 674800dd8cSBarry Smith 683a40ed3dSBarry Smith PetscFunctionBegin; 69fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 70fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 7139e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 72fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 73fbe28522SBarry Smith G = snes->work[1]; 746b5873e3SBarry Smith Ytmp = snes->work[2]; 754800dd8cSBarry Smith 764c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 774c49b128SBarry Smith snes->iter = 0; 784c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 797e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 804800dd8cSBarry Smith 815334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 82cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 830f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 84fbe28522SBarry Smith snes->norm = fnorm; 850f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 864800dd8cSBarry Smith delta = neP->delta0*fnorm; 874800dd8cSBarry Smith neP->delta = delta; 880462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 8994a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 904800dd8cSBarry Smith 91184914b5SBarry Smith if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 92b37302c6SLois Curfman McInnes 93d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 94d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 95d034289bSBarry Smith 96a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 9778b31e54SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 9878b31e54SBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 996831982aSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 100df60cc22SBarry Smith 1014800dd8cSBarry Smith for (i=0; i<maxits; i++) { 1025334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1035334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1042deea200SLois Curfman McInnes 1052deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10678b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr); 107329f5518SBarry Smith snes->linear_its += lits; 108af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 109cddf8d76SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr); 1102deea200SLois Curfman McInnes norm1 = norm; 1116b5873e3SBarry Smith while(1) { 112393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 1132deea200SLois Curfman McInnes norm = norm1; 1146b5873e3SBarry Smith 1152deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 116eafb4bcbSBarry Smith if (norm >= delta) { 117fbe28522SBarry Smith norm = delta/norm; 118fbe28522SBarry Smith gpnorm = (1.0 - norm)*fnorm; 119eafb4bcbSBarry Smith cnorm = norm; 120af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm); 121393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 122eafb4bcbSBarry Smith norm = gpnorm; 123fbe28522SBarry Smith ynorm = delta; 124fbe28522SBarry Smith } else { 125fbe28522SBarry Smith gpnorm = 0.0; 126af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n"); 127fbe28522SBarry Smith ynorm = norm; 128fbe28522SBarry Smith } 1292deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr); /* Y <- X - Y */ 13081b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1315334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 132cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1334800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 134fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1354800dd8cSBarry Smith 1364800dd8cSBarry Smith /* Update size of trust region */ 1374800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1384800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1394800dd8cSBarry Smith else delta *= neP->delta3; 140af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 141af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1424800dd8cSBarry Smith neP->delta = delta; 1434800dd8cSBarry Smith if (rho > neP->sigma) break; 144af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n"); 1456b5873e3SBarry Smith /* check to see if progress is hopeless */ 14652392280SLois Curfman McInnes neP->itflag = 0; 147184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 148184914b5SBarry Smith if (reason) { 14952392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 150454a90a3SBarry Smith breakout = 1; 151454a90a3SBarry Smith break; 15252392280SLois Curfman McInnes } 15352392280SLois Curfman McInnes snes->nfailures++; 1544800dd8cSBarry Smith } 1551acabf8cSLois Curfman McInnes if (!breakout) { 1564800dd8cSBarry Smith fnorm = gnorm; 1570f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 158c509a162SBarry Smith snes->iter = i+1; 159fbe28522SBarry Smith snes->norm = fnorm; 1600f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16139e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16239e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 163329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1640462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16594a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1664800dd8cSBarry Smith 1674800dd8cSBarry Smith /* Test for convergence */ 16852392280SLois Curfman McInnes neP->itflag = 1; 169184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 170184914b5SBarry Smith if (reason) { 17138442cffSBarry Smith break; 17238442cffSBarry Smith } 1731acabf8cSLois Curfman McInnes } else { 1741acabf8cSLois Curfman McInnes break; 1751acabf8cSLois Curfman McInnes } 17638442cffSBarry Smith } 17738442cffSBarry Smith /* Verify solution is in corect location */ 178e15f7bb5SBarry Smith if (X != snes->vec_sol) { 179393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 180e15f7bb5SBarry Smith } 181e15f7bb5SBarry Smith if (F != snes->vec_func) { 182e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 183e15f7bb5SBarry Smith } 18439e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18539e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18652392280SLois Curfman McInnes if (i == maxits) { 187af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits); 18852392280SLois Curfman McInnes i--; 189184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 19052392280SLois Curfman McInnes } 19152392280SLois Curfman McInnes *its = i+1; 1920f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 193184914b5SBarry Smith snes->reason = reason; 1940f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1953a40ed3dSBarry Smith PetscFunctionReturn(0); 1964800dd8cSBarry Smith } 1974800dd8cSBarry Smith /*------------------------------------------------------------*/ 1985615d1e5SSatish Balay #undef __FUNC__ 199b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_TR" 200f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_TR(SNES snes) 2014800dd8cSBarry Smith { 202fbe28522SBarry Smith int ierr; 2033a40ed3dSBarry Smith 2043a40ed3dSBarry Smith PetscFunctionBegin; 20581b6cf68SLois Curfman McInnes snes->nwork = 4; 206d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 207393d2d9aSLois Curfman McInnes PLogObjectParents(snes,snes->nwork,snes->work); 20881b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 2104800dd8cSBarry Smith } 2114800dd8cSBarry Smith /*------------------------------------------------------------*/ 2125615d1e5SSatish Balay #undef __FUNC__ 213b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_TR" 214e1311b90SBarry Smith static int SNESDestroy_EQ_TR(SNES snes) 2154800dd8cSBarry Smith { 216393d2d9aSLois Curfman McInnes int ierr; 2175baf8537SBarry Smith 2183a40ed3dSBarry Smith PetscFunctionBegin; 2195baf8537SBarry Smith if (snes->nwork) { 2204b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2215baf8537SBarry Smith } 222606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 2244800dd8cSBarry Smith } 2254800dd8cSBarry Smith /*------------------------------------------------------------*/ 2264800dd8cSBarry Smith 2275615d1e5SSatish Balay #undef __FUNC__ 228b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_TR" 229f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes) 2304800dd8cSBarry Smith { 2316831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 232f1af5d2fSBarry Smith int ierr; 2334800dd8cSBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 2354bbc92c1SBarry Smith ierr = OptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 236186905e3SBarry Smith ierr = OptionsDouble("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 237186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 238186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 239186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 240186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 241186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 242186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 243186905e3SBarry Smith ierr = OptionsDouble("-snes_eq_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 2444bbc92c1SBarry Smith ierr = OptionsTail();CHKERRQ(ierr); 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2464800dd8cSBarry Smith } 2474800dd8cSBarry Smith 2485615d1e5SSatish Balay #undef __FUNC__ 249b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_TR" 250e1311b90SBarry Smith static int SNESView_EQ_TR(SNES snes,Viewer viewer) 251a935fc98SLois Curfman McInnes { 2526831982aSBarry Smith SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data; 25351695f54SLois Curfman McInnes int ierr; 2546831982aSBarry Smith PetscTruth isascii; 255a935fc98SLois Curfman McInnes 2563a40ed3dSBarry Smith PetscFunctionBegin; 2576831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 2580f5bd95cSBarry Smith if (isascii) { 259d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 260d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2615cd90555SBarry Smith } else { 262*29bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 26319bcc07fSBarry Smith } 2643a40ed3dSBarry Smith PetscFunctionReturn(0); 265a935fc98SLois Curfman McInnes } 266a935fc98SLois Curfman McInnes 26752392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2685615d1e5SSatish Balay #undef __FUNC__ 269b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESConverged_EQ_TR" 2708ac28fe4SSatish Balay /*@C 271f525115eSLois Curfman McInnes SNESConverged_EQ_TR - Monitors the convergence of the trust region 2726831982aSBarry Smith method SNESEQTR for solving systems of nonlinear equations (default). 27352392280SLois Curfman McInnes 274c7afd0dbSLois Curfman McInnes Collective on SNES 275c7afd0dbSLois Curfman McInnes 27652392280SLois Curfman McInnes Input Parameters: 277c7afd0dbSLois Curfman McInnes + snes - the SNES context 27852392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 27952392280SLois Curfman McInnes . pnorm - 2-norm of current step 28052392280SLois Curfman McInnes . fnorm - 2-norm of function 281c7afd0dbSLois Curfman McInnes - dummy - unused context 282fee21e36SBarry Smith 283184914b5SBarry Smith Output Parameter: 284184914b5SBarry Smith . reason - one of 285184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol), 286184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 287184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 288184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 289184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 290184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 291184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 29252392280SLois Curfman McInnes 29352392280SLois Curfman McInnes where 294c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 295c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 296c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 297c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 298c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 299c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 300c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 301c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 302c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 303c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 30452392280SLois Curfman McInnes 30515091d37SBarry Smith Level: intermediate 30615091d37SBarry Smith 30752392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 30852392280SLois Curfman McInnes 30952392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 31052392280SLois Curfman McInnes @*/ 311184914b5SBarry Smith int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy) 31252392280SLois Curfman McInnes { 3136831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data; 314184914b5SBarry Smith int ierr; 315ddd12767SBarry Smith 3163a40ed3dSBarry Smith PetscFunctionBegin; 3173a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 318*29bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 3193a40ed3dSBarry Smith } 32052392280SLois Curfman McInnes 321d252947aSBarry Smith if (fnorm != fnorm) { 322981c4779SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n"); 323184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 324184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 325184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 326184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 327184914b5SBarry Smith } else if (neP->itflag) { 328184914b5SBarry Smith ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3291acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 330184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs); 331184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 332184914b5SBarry Smith } else { 333184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 33452392280SLois Curfman McInnes } 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 33652392280SLois Curfman McInnes } 33752392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 338fb2e594dSBarry Smith EXTERN_C_BEGIN 3395615d1e5SSatish Balay #undef __FUNC__ 340b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_TR" 341f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes) 3424800dd8cSBarry Smith { 3436831982aSBarry Smith SNES_EQ_TR *neP; 3444800dd8cSBarry Smith 3453a40ed3dSBarry Smith PetscFunctionBegin; 3463a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 347*29bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 3483a40ed3dSBarry Smith } 349f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_TR; 350f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_TR; 351f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_TR; 352f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_TR; 353f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_TR; 354f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_TR; 3555baf8537SBarry Smith snes->nwork = 0; 356fbe28522SBarry Smith 3576831982aSBarry Smith neP = PetscNew(SNES_EQ_TR);CHKPTRQ(neP); 3586831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_TR)); 359fbe28522SBarry Smith snes->data = (void*)neP; 360fbe28522SBarry Smith neP->mu = 0.25; 361fbe28522SBarry Smith neP->eta = 0.75; 362fbe28522SBarry Smith neP->delta = 0.0; 363fbe28522SBarry Smith neP->delta0 = 0.2; 364fbe28522SBarry Smith neP->delta1 = 0.3; 365fbe28522SBarry Smith neP->delta2 = 0.75; 366fbe28522SBarry Smith neP->delta3 = 2.0; 367fbe28522SBarry Smith neP->sigma = 0.0001; 368fbe28522SBarry Smith neP->itflag = 0; 36952392280SLois Curfman McInnes neP->rnorm0 = 0; 37052392280SLois Curfman McInnes neP->ttol = 0; 3713a40ed3dSBarry Smith PetscFunctionReturn(0); 3724800dd8cSBarry Smith } 373fb2e594dSBarry Smith EXTERN_C_END 37482bf6240SBarry Smith 375