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__ 104b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private" 114b27c08aSLois Curfman McInnes int SNES_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; 154b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_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) { 284b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNES_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) { 354b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 364b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNES_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 /* 434b27c08aSLois Curfman McInnes SNESSolve_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__ 554b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR" 564b27c08aSLois Curfman McInnes static int SNESSolve_TR(SNES snes,int *its) 574800dd8cSBarry Smith { 584b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_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; 65184914b5SBarry Smith SNESConvergedReason reason; 667c4af3dcSLois Curfman McInnes PetscTruth conv; 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); 90*94b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 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) { 1004b27c08aSLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 1014b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Using Krylov convergence test SNES_TR_KSPConverged_Private\n"); 1027c4af3dcSLois Curfman McInnes } 103df60cc22SBarry Smith 1044800dd8cSBarry Smith for (i=0; i<maxits; i++) { 1055334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 106*94b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1072deea200SLois Curfman McInnes 1082deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 109*94b7f48cSBarry Smith ierr = KSPSetRhs(snes->ksp,F);CHKERRQ(ierr); 110*94b7f48cSBarry Smith ierr = KSPSetSolution(snes->ksp,Ytmp);CHKERRQ(ierr); 111*94b7f48cSBarry Smith ierr = KSPSolve(snes->ksp);CHKERRQ(ierr); 112c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 113329f5518SBarry Smith snes->linear_its += lits; 1144b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 115064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 116064f8208SBarry Smith norm1 = nrm; 1176b5873e3SBarry Smith while(1) { 118393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 119064f8208SBarry Smith nrm = norm1; 1206b5873e3SBarry Smith 1212deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 122064f8208SBarry Smith if (nrm >= delta) { 123064f8208SBarry Smith nrm = delta/nrm; 124064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 125064f8208SBarry Smith cnorm = nrm; 1264b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Scaling direction by %g\n",nrm); 127393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 128064f8208SBarry Smith nrm = gpnorm; 129fbe28522SBarry Smith ynorm = delta; 130fbe28522SBarry Smith } else { 131fbe28522SBarry Smith gpnorm = 0.0; 1324b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Direction is in Trust Region\n"); 133064f8208SBarry Smith ynorm = nrm; 134fbe28522SBarry Smith } 1352deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr); /* Y <- X - Y */ 13681b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1375334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 138cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1394800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 140fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1414800dd8cSBarry Smith 1424800dd8cSBarry Smith /* Update size of trust region */ 1434800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1444800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1454800dd8cSBarry Smith else delta *= neP->delta3; 1464b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 1474b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1484800dd8cSBarry Smith neP->delta = delta; 1494800dd8cSBarry Smith if (rho > neP->sigma) break; 1504b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Trying again in smaller region\n"); 1516b5873e3SBarry Smith /* check to see if progress is hopeless */ 15252392280SLois Curfman McInnes neP->itflag = 0; 153184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 154184914b5SBarry Smith if (reason) { 15552392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 1565ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 157454a90a3SBarry Smith breakout = 1; 158454a90a3SBarry Smith break; 15952392280SLois Curfman McInnes } 16050ffb88aSMatthew Knepley snes->numFailures++; 1614800dd8cSBarry Smith } 1621acabf8cSLois Curfman McInnes if (!breakout) { 1634800dd8cSBarry Smith fnorm = gnorm; 1640f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 165c509a162SBarry Smith snes->iter = i+1; 166fbe28522SBarry Smith snes->norm = fnorm; 1670f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16839e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16939e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 170329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1710462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 17294a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1734800dd8cSBarry Smith 1744800dd8cSBarry Smith /* Test for convergence */ 17552392280SLois Curfman McInnes neP->itflag = 1; 176184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 177184914b5SBarry Smith if (reason) { 17838442cffSBarry Smith break; 17938442cffSBarry Smith } 1801acabf8cSLois Curfman McInnes } else { 1811acabf8cSLois Curfman McInnes break; 1821acabf8cSLois Curfman McInnes } 18338442cffSBarry Smith } 18438442cffSBarry Smith /* Verify solution is in corect location */ 185e15f7bb5SBarry Smith if (X != snes->vec_sol) { 186393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 187e15f7bb5SBarry Smith } 188e15f7bb5SBarry Smith if (F != snes->vec_func) { 189e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 190e15f7bb5SBarry Smith } 19139e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 19239e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 19352392280SLois Curfman McInnes if (i == maxits) { 1944b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Maximum number of iterations has been reached: %d\n",maxits); 19552392280SLois Curfman McInnes i--; 196184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 19752392280SLois Curfman McInnes } 19852392280SLois Curfman McInnes *its = i+1; 1990f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 200184914b5SBarry Smith snes->reason = reason; 2010f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 2034800dd8cSBarry Smith } 2044800dd8cSBarry Smith /*------------------------------------------------------------*/ 2054a2ae208SSatish Balay #undef __FUNCT__ 2064b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR" 2074b27c08aSLois Curfman McInnes static int SNESSetUp_TR(SNES snes) 2084800dd8cSBarry Smith { 209fbe28522SBarry Smith int ierr; 2103a40ed3dSBarry Smith 2113a40ed3dSBarry Smith PetscFunctionBegin; 21281b6cf68SLois Curfman McInnes snes->nwork = 4; 213d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 214b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 21581b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2163a40ed3dSBarry Smith PetscFunctionReturn(0); 2174800dd8cSBarry Smith } 2184800dd8cSBarry Smith /*------------------------------------------------------------*/ 2194a2ae208SSatish Balay #undef __FUNCT__ 2204b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR" 2214b27c08aSLois Curfman McInnes static int SNESDestroy_TR(SNES snes) 2224800dd8cSBarry Smith { 223393d2d9aSLois Curfman McInnes int ierr; 2245baf8537SBarry Smith 2253a40ed3dSBarry Smith PetscFunctionBegin; 2265baf8537SBarry Smith if (snes->nwork) { 2274b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2285baf8537SBarry Smith } 229606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2303a40ed3dSBarry Smith PetscFunctionReturn(0); 2314800dd8cSBarry Smith } 2324800dd8cSBarry Smith /*------------------------------------------------------------*/ 2334800dd8cSBarry Smith 2344a2ae208SSatish Balay #undef __FUNCT__ 2354b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR" 2364b27c08aSLois Curfman McInnes static int SNESSetFromOptions_TR(SNES snes) 2374800dd8cSBarry Smith { 2384b27c08aSLois Curfman McInnes SNES_TR *ctx = (SNES_TR *)snes->data; 239f1af5d2fSBarry Smith int ierr; 2404800dd8cSBarry Smith 2413a40ed3dSBarry Smith PetscFunctionBegin; 242b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 24387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 2444b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 2454b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 2464b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 2474b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 2484b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 2494b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 2504b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 251b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2523a40ed3dSBarry Smith PetscFunctionReturn(0); 2534800dd8cSBarry Smith } 2544800dd8cSBarry Smith 2554a2ae208SSatish Balay #undef __FUNCT__ 2564b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR" 2574b27c08aSLois Curfman McInnes static int SNESView_TR(SNES snes,PetscViewer viewer) 258a935fc98SLois Curfman McInnes { 2594b27c08aSLois Curfman McInnes SNES_TR *tr = (SNES_TR *)snes->data; 26051695f54SLois Curfman McInnes int ierr; 2616831982aSBarry Smith PetscTruth isascii; 262a935fc98SLois Curfman McInnes 2633a40ed3dSBarry Smith PetscFunctionBegin; 264b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 2650f5bd95cSBarry Smith if (isascii) { 266b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 267b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2685cd90555SBarry Smith } else { 26929bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 27019bcc07fSBarry Smith } 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 272a935fc98SLois Curfman McInnes } 273a935fc98SLois Curfman McInnes 27452392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2754a2ae208SSatish Balay #undef __FUNCT__ 2764b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR" 2778ac28fe4SSatish Balay /*@C 2784b27c08aSLois Curfman McInnes SNESConverged_TR - Monitors the convergence of the trust region 2794b27c08aSLois Curfman McInnes method SNESTR for solving systems of nonlinear equations (default). 28052392280SLois Curfman McInnes 281c7afd0dbSLois Curfman McInnes Collective on SNES 282c7afd0dbSLois Curfman McInnes 28352392280SLois Curfman McInnes Input Parameters: 284c7afd0dbSLois Curfman McInnes + snes - the SNES context 28552392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 28652392280SLois Curfman McInnes . pnorm - 2-norm of current step 28752392280SLois Curfman McInnes . fnorm - 2-norm of function 288c7afd0dbSLois Curfman McInnes - dummy - unused context 289fee21e36SBarry Smith 290184914b5SBarry Smith Output Parameter: 291184914b5SBarry Smith . reason - one of 292184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol), 293184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 294184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 295184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 296184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 297184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 298184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 29952392280SLois Curfman McInnes 30052392280SLois Curfman McInnes where 301c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 302c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 303c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 304c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 305c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 306c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 307c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 308c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 309c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 310c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 31152392280SLois Curfman McInnes 31215091d37SBarry Smith Level: intermediate 31315091d37SBarry Smith 31452392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 31552392280SLois Curfman McInnes 31652392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 31752392280SLois Curfman McInnes @*/ 3184b27c08aSLois Curfman McInnes int SNESConverged_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 31952392280SLois Curfman McInnes { 3204b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR *)snes->data; 321184914b5SBarry Smith int ierr; 322ddd12767SBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 324d252947aSBarry Smith if (fnorm != fnorm) { 3254b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR:Failed to converged, function norm is NaN\n"); 326184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 327184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 3284b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 329184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 330184914b5SBarry Smith } else if (neP->itflag) { 3314b27c08aSLois Curfman McInnes ierr = SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3321acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 3334b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs); 334184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 335184914b5SBarry Smith } else { 336184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 33752392280SLois Curfman McInnes } 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 33952392280SLois Curfman McInnes } 34052392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 341fb2e594dSBarry Smith EXTERN_C_BEGIN 3424a2ae208SSatish Balay #undef __FUNCT__ 3434b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR" 3444b27c08aSLois Curfman McInnes int SNESCreate_TR(SNES snes) 3454800dd8cSBarry Smith { 3464b27c08aSLois Curfman McInnes SNES_TR *neP; 34782502324SSatish Balay int ierr; 3484800dd8cSBarry Smith 3493a40ed3dSBarry Smith PetscFunctionBegin; 3504b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_TR; 3514b27c08aSLois Curfman McInnes snes->solve = SNESSolve_TR; 3524b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_TR; 3534b27c08aSLois Curfman McInnes snes->converged = SNESConverged_TR; 3544b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_TR; 3554b27c08aSLois Curfman McInnes snes->view = SNESView_TR; 3565baf8537SBarry Smith snes->nwork = 0; 357fbe28522SBarry Smith 3584b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_TR,&neP);CHKERRQ(ierr); 3594b27c08aSLois Curfman McInnes PetscLogObjectMemory(snes,sizeof(SNES_TR)); 360fbe28522SBarry Smith snes->data = (void*)neP; 361fbe28522SBarry Smith neP->mu = 0.25; 362fbe28522SBarry Smith neP->eta = 0.75; 363fbe28522SBarry Smith neP->delta = 0.0; 364fbe28522SBarry Smith neP->delta0 = 0.2; 365fbe28522SBarry Smith neP->delta1 = 0.3; 366fbe28522SBarry Smith neP->delta2 = 0.75; 367fbe28522SBarry Smith neP->delta3 = 2.0; 368fbe28522SBarry Smith neP->sigma = 0.0001; 369fbe28522SBarry Smith neP->itflag = 0; 37052392280SLois Curfman McInnes neP->rnorm0 = 0; 37152392280SLois Curfman McInnes neP->ttol = 0; 3723a40ed3dSBarry Smith PetscFunctionReturn(0); 3734800dd8cSBarry Smith } 374fb2e594dSBarry Smith EXTERN_C_END 37582bf6240SBarry Smith 376