173f4d377SMatthew Knepley /*$Id: tr.c,v 1.128 2001/08/07 03:04:11 balay Exp $*/ 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__ 104a2ae208SSatish Balay #define __FUNCT__ "SNES_EQ_TR_KSPConverged_Private" 1187828ca2SBarry Smith int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,PetscReal 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; 17064f8208SBarry Smith PetscReal nrm; 183acb151aSSatish Balay int ierr; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 2198120f82SLois Curfman McInnes if (snes->ksp_ewconv) { 2229bbc08cSBarry 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) { 28b0a32e0cSBarry Smith PetscLogInfo(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); 33064f8208SBarry Smith ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 34064f8208SBarry Smith if (nrm >= neP->delta) { 35b0a32e0cSBarry Smith PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 36064f8208SBarry Smith PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,nrm); 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 */ 544a2ae208SSatish Balay #undef __FUNCT__ 554a2ae208SSatish Balay #define __FUNCT__ "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; 62064f8208SBarry Smith PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1; 63ea709b57SSatish Balay PetscScalar mone = -1.0,cnorm; 64df60cc22SBarry Smith KSP ksp; 65df60cc22SBarry Smith SLES sles; 66184914b5SBarry Smith SNESConvergedReason reason; 677c4af3dcSLois Curfman McInnes PetscTruth conv; 684800dd8cSBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 70fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 71fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 7239e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 73fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 74fbe28522SBarry Smith G = snes->work[1]; 756b5873e3SBarry Smith Ytmp = snes->work[2]; 764800dd8cSBarry Smith 774c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 784c49b128SBarry Smith snes->iter = 0; 794c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 807e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 814800dd8cSBarry Smith 825334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 83cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 840f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 85fbe28522SBarry Smith snes->norm = fnorm; 860f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 874800dd8cSBarry Smith delta = neP->delta0*fnorm; 884800dd8cSBarry Smith neP->delta = delta; 890462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 9094a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 914800dd8cSBarry Smith 92184914b5SBarry Smith if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 93b37302c6SLois Curfman McInnes 94d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 95d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 96d034289bSBarry Smith 97a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 987c4af3dcSLois Curfman McInnes ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 997c4af3dcSLois Curfman McInnes if (!conv) { 10078b31e54SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 10178b31e54SBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 1026831982aSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 103a85fdd5bSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_EQ_TR: Using Krylov convergence test SNES_EQ_TR_KSPConverged_Private\n"); 1047c4af3dcSLois Curfman McInnes } 105df60cc22SBarry Smith 1064800dd8cSBarry Smith for (i=0; i<maxits; i++) { 1075334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1085334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1092deea200SLois Curfman McInnes 1102deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 11178b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr); 112329f5518SBarry Smith snes->linear_its += lits; 113b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 114064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 115064f8208SBarry Smith norm1 = nrm; 1166b5873e3SBarry Smith while(1) { 117393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 118064f8208SBarry Smith nrm = norm1; 1196b5873e3SBarry Smith 1202deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 121064f8208SBarry Smith if (nrm >= delta) { 122064f8208SBarry Smith nrm = delta/nrm; 123064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 124064f8208SBarry Smith cnorm = nrm; 125064f8208SBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",nrm); 126393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 127064f8208SBarry Smith nrm = gpnorm; 128fbe28522SBarry Smith ynorm = delta; 129fbe28522SBarry Smith } else { 130fbe28522SBarry Smith gpnorm = 0.0; 131b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n"); 132064f8208SBarry Smith ynorm = nrm; 133fbe28522SBarry Smith } 1342deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr); /* Y <- X - Y */ 13581b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1365334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 137cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1384800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 139fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1404800dd8cSBarry Smith 1414800dd8cSBarry Smith /* Update size of trust region */ 1424800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1434800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1444800dd8cSBarry Smith else delta *= neP->delta3; 145b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 146b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1474800dd8cSBarry Smith neP->delta = delta; 1484800dd8cSBarry Smith if (rho > neP->sigma) break; 149b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n"); 1506b5873e3SBarry Smith /* check to see if progress is hopeless */ 15152392280SLois Curfman McInnes neP->itflag = 0; 152184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 153184914b5SBarry Smith if (reason) { 15452392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 155*5ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 156454a90a3SBarry Smith breakout = 1; 157454a90a3SBarry Smith break; 15852392280SLois Curfman McInnes } 15950ffb88aSMatthew Knepley snes->numFailures++; 1604800dd8cSBarry Smith } 1611acabf8cSLois Curfman McInnes if (!breakout) { 1624800dd8cSBarry Smith fnorm = gnorm; 1630f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 164c509a162SBarry Smith snes->iter = i+1; 165fbe28522SBarry Smith snes->norm = fnorm; 1660f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16739e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16839e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 169329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1700462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 17194a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1724800dd8cSBarry Smith 1734800dd8cSBarry Smith /* Test for convergence */ 17452392280SLois Curfman McInnes neP->itflag = 1; 175184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 176184914b5SBarry Smith if (reason) { 17738442cffSBarry Smith break; 17838442cffSBarry Smith } 1791acabf8cSLois Curfman McInnes } else { 1801acabf8cSLois Curfman McInnes break; 1811acabf8cSLois Curfman McInnes } 18238442cffSBarry Smith } 18338442cffSBarry Smith /* Verify solution is in corect location */ 184e15f7bb5SBarry Smith if (X != snes->vec_sol) { 185393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 186e15f7bb5SBarry Smith } 187e15f7bb5SBarry Smith if (F != snes->vec_func) { 188e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 189e15f7bb5SBarry Smith } 19039e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 19139e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 19252392280SLois Curfman McInnes if (i == maxits) { 193b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits); 19452392280SLois Curfman McInnes i--; 195184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 19652392280SLois Curfman McInnes } 19752392280SLois Curfman McInnes *its = i+1; 1980f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 199184914b5SBarry Smith snes->reason = reason; 2000f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2013a40ed3dSBarry Smith PetscFunctionReturn(0); 2024800dd8cSBarry Smith } 2034800dd8cSBarry Smith /*------------------------------------------------------------*/ 2044a2ae208SSatish Balay #undef __FUNCT__ 2054a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_TR" 206f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_TR(SNES snes) 2074800dd8cSBarry Smith { 208fbe28522SBarry Smith int ierr; 2093a40ed3dSBarry Smith 2103a40ed3dSBarry Smith PetscFunctionBegin; 21181b6cf68SLois Curfman McInnes snes->nwork = 4; 212d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 213b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 21481b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 2164800dd8cSBarry Smith } 2174800dd8cSBarry Smith /*------------------------------------------------------------*/ 2184a2ae208SSatish Balay #undef __FUNCT__ 2194a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_TR" 220e1311b90SBarry Smith static int SNESDestroy_EQ_TR(SNES snes) 2214800dd8cSBarry Smith { 222393d2d9aSLois Curfman McInnes int ierr; 2235baf8537SBarry Smith 2243a40ed3dSBarry Smith PetscFunctionBegin; 2255baf8537SBarry Smith if (snes->nwork) { 2264b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2275baf8537SBarry Smith } 228606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 2304800dd8cSBarry Smith } 2314800dd8cSBarry Smith /*------------------------------------------------------------*/ 2324800dd8cSBarry Smith 2334a2ae208SSatish Balay #undef __FUNCT__ 2344a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_TR" 235f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes) 2364800dd8cSBarry Smith { 2376831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 238f1af5d2fSBarry Smith int ierr; 2394800dd8cSBarry Smith 2403a40ed3dSBarry Smith PetscFunctionBegin; 241b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 24287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 24387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 24487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 24587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 24687828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 24787828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 24887828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 24987828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 250b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2513a40ed3dSBarry Smith PetscFunctionReturn(0); 2524800dd8cSBarry Smith } 2534800dd8cSBarry Smith 2544a2ae208SSatish Balay #undef __FUNCT__ 2554a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_TR" 256b0a32e0cSBarry Smith static int SNESView_EQ_TR(SNES snes,PetscViewer viewer) 257a935fc98SLois Curfman McInnes { 2586831982aSBarry Smith SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data; 25951695f54SLois Curfman McInnes int ierr; 2606831982aSBarry Smith PetscTruth isascii; 261a935fc98SLois Curfman McInnes 2623a40ed3dSBarry Smith PetscFunctionBegin; 263b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 2640f5bd95cSBarry Smith if (isascii) { 265b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 266b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2675cd90555SBarry Smith } else { 26829bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 26919bcc07fSBarry Smith } 2703a40ed3dSBarry Smith PetscFunctionReturn(0); 271a935fc98SLois Curfman McInnes } 272a935fc98SLois Curfman McInnes 27352392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "SNESConverged_EQ_TR" 2768ac28fe4SSatish Balay /*@C 277f525115eSLois Curfman McInnes SNESConverged_EQ_TR - Monitors the convergence of the trust region 2786831982aSBarry Smith method SNESEQTR for solving systems of nonlinear equations (default). 27952392280SLois Curfman McInnes 280c7afd0dbSLois Curfman McInnes Collective on SNES 281c7afd0dbSLois Curfman McInnes 28252392280SLois Curfman McInnes Input Parameters: 283c7afd0dbSLois Curfman McInnes + snes - the SNES context 28452392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 28552392280SLois Curfman McInnes . pnorm - 2-norm of current step 28652392280SLois Curfman McInnes . fnorm - 2-norm of function 287c7afd0dbSLois Curfman McInnes - dummy - unused context 288fee21e36SBarry Smith 289184914b5SBarry Smith Output Parameter: 290184914b5SBarry Smith . reason - one of 291184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol), 292184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 293184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 294184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 295184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 296184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 297184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 29852392280SLois Curfman McInnes 29952392280SLois Curfman McInnes where 300c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 301c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 302c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 303c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 304c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 305c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 306c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 307c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 308c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 309c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 31052392280SLois Curfman McInnes 31115091d37SBarry Smith Level: intermediate 31215091d37SBarry Smith 31352392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 31452392280SLois Curfman McInnes 31552392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 31652392280SLois Curfman McInnes @*/ 31787828ca2SBarry Smith int SNESConverged_EQ_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 31852392280SLois Curfman McInnes { 3196831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data; 320184914b5SBarry Smith int ierr; 321ddd12767SBarry Smith 3223a40ed3dSBarry Smith PetscFunctionBegin; 3233a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 32429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 3253a40ed3dSBarry Smith } 32652392280SLois Curfman McInnes 327d252947aSBarry Smith if (fnorm != fnorm) { 328b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n"); 329184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 330184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 331b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 332184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 333184914b5SBarry Smith } else if (neP->itflag) { 334184914b5SBarry Smith ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3351acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 336b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs); 337184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 338184914b5SBarry Smith } else { 339184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 34052392280SLois Curfman McInnes } 3413a40ed3dSBarry Smith PetscFunctionReturn(0); 34252392280SLois Curfman McInnes } 34352392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 344fb2e594dSBarry Smith EXTERN_C_BEGIN 3454a2ae208SSatish Balay #undef __FUNCT__ 3464a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_TR" 347f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes) 3484800dd8cSBarry Smith { 3496831982aSBarry Smith SNES_EQ_TR *neP; 35082502324SSatish Balay int ierr; 3514800dd8cSBarry Smith 3523a40ed3dSBarry Smith PetscFunctionBegin; 3533a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 35429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 3553a40ed3dSBarry Smith } 356f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_TR; 357f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_TR; 358f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_TR; 359f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_TR; 360f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_TR; 361f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_TR; 3625baf8537SBarry Smith snes->nwork = 0; 363fbe28522SBarry Smith 364b0a32e0cSBarry Smith ierr = PetscNew(SNES_EQ_TR,&neP);CHKERRQ(ierr); 365b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(SNES_EQ_TR)); 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; 375fbe28522SBarry Smith neP->itflag = 0; 37652392280SLois Curfman McInnes neP->rnorm0 = 0; 37752392280SLois Curfman McInnes neP->ttol = 0; 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 3794800dd8cSBarry Smith } 380fb2e594dSBarry Smith EXTERN_C_END 38182bf6240SBarry Smith 382