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 46ddbbdb52SLois Curfman McInnes 474800dd8cSBarry Smith */ 484a2ae208SSatish Balay #undef __FUNCT__ 494b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR" 50c9780b6fSBarry Smith static int SNESSolve_TR(SNES snes) 514800dd8cSBarry Smith { 524b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR*)snes->data; 536b5873e3SBarry Smith Vec X,F,Y,G,TMP,Ytmp; 540462333dSBarry Smith int maxits,i,ierr,lits,breakout = 0; 55112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 56064f8208SBarry Smith PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1; 57ea709b57SSatish Balay PetscScalar mone = -1.0,cnorm; 58df60cc22SBarry Smith KSP ksp; 59184914b5SBarry Smith SNESConvergedReason reason; 607c4af3dcSLois Curfman McInnes PetscTruth conv; 614800dd8cSBarry Smith 623a40ed3dSBarry Smith PetscFunctionBegin; 63fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 64fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 6539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 66fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 67fbe28522SBarry Smith G = snes->work[1]; 686b5873e3SBarry Smith Ytmp = snes->work[2]; 694800dd8cSBarry Smith 704c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 714c49b128SBarry Smith snes->iter = 0; 724c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 737e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 744800dd8cSBarry Smith 755334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 76cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 770f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 78fbe28522SBarry Smith snes->norm = fnorm; 790f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 804800dd8cSBarry Smith delta = neP->delta0*fnorm; 814800dd8cSBarry Smith neP->delta = delta; 820462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 8394a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 8494b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 854800dd8cSBarry Smith 86c9780b6fSBarry Smith if (fnorm < snes->atol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 87b37302c6SLois Curfman McInnes 88d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 89d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 90d034289bSBarry Smith 91a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 927c4af3dcSLois Curfman McInnes ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 937c4af3dcSLois Curfman McInnes if (!conv) { 944b27c08aSLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 954b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Using Krylov convergence test SNES_TR_KSPConverged_Private\n"); 967c4af3dcSLois Curfman McInnes } 97df60cc22SBarry Smith 984800dd8cSBarry Smith for (i=0; i<maxits; i++) { 995334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 10094b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1012deea200SLois Curfman McInnes 1022deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10394b7f48cSBarry Smith ierr = KSPSetRhs(snes->ksp,F);CHKERRQ(ierr); 10494b7f48cSBarry Smith ierr = KSPSetSolution(snes->ksp,Ytmp);CHKERRQ(ierr); 10594b7f48cSBarry Smith ierr = KSPSolve(snes->ksp);CHKERRQ(ierr); 106c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 107329f5518SBarry Smith snes->linear_its += lits; 1084b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 109064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 110064f8208SBarry Smith norm1 = nrm; 1116b5873e3SBarry Smith while(1) { 112393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 113064f8208SBarry Smith nrm = norm1; 1146b5873e3SBarry Smith 1152deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 116064f8208SBarry Smith if (nrm >= delta) { 117064f8208SBarry Smith nrm = delta/nrm; 118064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 119064f8208SBarry Smith cnorm = nrm; 1204b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Scaling direction by %g\n",nrm); 121393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 122064f8208SBarry Smith nrm = gpnorm; 123fbe28522SBarry Smith ynorm = delta; 124fbe28522SBarry Smith } else { 125fbe28522SBarry Smith gpnorm = 0.0; 1264b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Direction is in Trust Region\n"); 127064f8208SBarry Smith ynorm = nrm; 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; 1404b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 1414b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1424800dd8cSBarry Smith neP->delta = delta; 1434800dd8cSBarry Smith if (rho > neP->sigma) break; 1444b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_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 */ 1505ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 151454a90a3SBarry Smith breakout = 1; 152454a90a3SBarry Smith break; 15352392280SLois Curfman McInnes } 15450ffb88aSMatthew Knepley snes->numFailures++; 1554800dd8cSBarry Smith } 1561acabf8cSLois Curfman McInnes if (!breakout) { 1574800dd8cSBarry Smith fnorm = gnorm; 1580f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 159c509a162SBarry Smith snes->iter = i+1; 160fbe28522SBarry Smith snes->norm = fnorm; 1610f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16239e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16339e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 164329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1650462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16694a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1674800dd8cSBarry Smith 1684800dd8cSBarry Smith /* Test for convergence */ 16952392280SLois Curfman McInnes neP->itflag = 1; 170184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 171184914b5SBarry Smith if (reason) { 17238442cffSBarry Smith break; 17338442cffSBarry Smith } 1741acabf8cSLois Curfman McInnes } else { 1751acabf8cSLois Curfman McInnes break; 1761acabf8cSLois Curfman McInnes } 17738442cffSBarry Smith } 17838442cffSBarry Smith /* Verify solution is in corect location */ 179e15f7bb5SBarry Smith if (X != snes->vec_sol) { 180393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 181e15f7bb5SBarry Smith } 182e15f7bb5SBarry Smith if (F != snes->vec_func) { 183e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 184e15f7bb5SBarry Smith } 18539e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18639e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18752392280SLois Curfman McInnes if (i == maxits) { 1884b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_TR: Maximum number of iterations has been reached: %d\n",maxits); 189184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 19052392280SLois Curfman McInnes } 1910f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 192184914b5SBarry Smith snes->reason = reason; 1930f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1943a40ed3dSBarry Smith PetscFunctionReturn(0); 1954800dd8cSBarry Smith } 1964800dd8cSBarry Smith /*------------------------------------------------------------*/ 1974a2ae208SSatish Balay #undef __FUNCT__ 1984b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR" 1994b27c08aSLois Curfman McInnes static int SNESSetUp_TR(SNES snes) 2004800dd8cSBarry Smith { 201fbe28522SBarry Smith int ierr; 2023a40ed3dSBarry Smith 2033a40ed3dSBarry Smith PetscFunctionBegin; 20481b6cf68SLois Curfman McInnes snes->nwork = 4; 205d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 206b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 20781b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 2094800dd8cSBarry Smith } 2104800dd8cSBarry Smith /*------------------------------------------------------------*/ 2114a2ae208SSatish Balay #undef __FUNCT__ 2124b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR" 2134b27c08aSLois Curfman McInnes static int SNESDestroy_TR(SNES snes) 2144800dd8cSBarry Smith { 215393d2d9aSLois Curfman McInnes int ierr; 2165baf8537SBarry Smith 2173a40ed3dSBarry Smith PetscFunctionBegin; 2185baf8537SBarry Smith if (snes->nwork) { 2194b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2205baf8537SBarry Smith } 221606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2223a40ed3dSBarry Smith PetscFunctionReturn(0); 2234800dd8cSBarry Smith } 2244800dd8cSBarry Smith /*------------------------------------------------------------*/ 2254800dd8cSBarry Smith 2264a2ae208SSatish Balay #undef __FUNCT__ 2274b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR" 2284b27c08aSLois Curfman McInnes static int SNESSetFromOptions_TR(SNES snes) 2294800dd8cSBarry Smith { 2304b27c08aSLois Curfman McInnes SNES_TR *ctx = (SNES_TR *)snes->data; 231f1af5d2fSBarry Smith int ierr; 2324800dd8cSBarry Smith 2333a40ed3dSBarry Smith PetscFunctionBegin; 234b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 23587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 2364b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 2374b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 2384b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 2394b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 2404b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 2414b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 2424b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 243b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2454800dd8cSBarry Smith } 2464800dd8cSBarry Smith 2474a2ae208SSatish Balay #undef __FUNCT__ 2484b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR" 2494b27c08aSLois Curfman McInnes static int SNESView_TR(SNES snes,PetscViewer viewer) 250a935fc98SLois Curfman McInnes { 2514b27c08aSLois Curfman McInnes SNES_TR *tr = (SNES_TR *)snes->data; 25251695f54SLois Curfman McInnes int ierr; 2536831982aSBarry Smith PetscTruth isascii; 254a935fc98SLois Curfman McInnes 2553a40ed3dSBarry Smith PetscFunctionBegin; 256b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 2570f5bd95cSBarry Smith if (isascii) { 258b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 259b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2605cd90555SBarry Smith } else { 26129bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 26219bcc07fSBarry Smith } 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 264a935fc98SLois Curfman McInnes } 265a935fc98SLois Curfman McInnes 26652392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2674a2ae208SSatish Balay #undef __FUNCT__ 2684b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR" 2698ac28fe4SSatish Balay /*@C 2704b27c08aSLois Curfman McInnes SNESConverged_TR - Monitors the convergence of the trust region 2714b27c08aSLois Curfman McInnes method SNESTR for solving systems of nonlinear equations (default). 27252392280SLois Curfman McInnes 273c7afd0dbSLois Curfman McInnes Collective on SNES 274c7afd0dbSLois Curfman McInnes 27552392280SLois Curfman McInnes Input Parameters: 276c7afd0dbSLois Curfman McInnes + snes - the SNES context 27752392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 27852392280SLois Curfman McInnes . pnorm - 2-norm of current step 27952392280SLois Curfman McInnes . fnorm - 2-norm of function 280c7afd0dbSLois Curfman McInnes - dummy - unused context 281fee21e36SBarry Smith 282184914b5SBarry Smith Output Parameter: 283184914b5SBarry Smith . reason - one of 284184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol), 285184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 286184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 287184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 288184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 289184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 290184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 29152392280SLois Curfman McInnes 29252392280SLois Curfman McInnes where 293c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 294c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 295c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 296c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 297c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 298c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 299c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 300c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 301c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 302c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 30352392280SLois Curfman McInnes 30415091d37SBarry Smith Level: intermediate 30515091d37SBarry Smith 30652392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 30752392280SLois Curfman McInnes 30852392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 30952392280SLois Curfman McInnes @*/ 3104b27c08aSLois Curfman McInnes int SNESConverged_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 31152392280SLois Curfman McInnes { 3124b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR *)snes->data; 313184914b5SBarry Smith int ierr; 314ddd12767SBarry Smith 3153a40ed3dSBarry Smith PetscFunctionBegin; 316d252947aSBarry Smith if (fnorm != fnorm) { 3174b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR:Failed to converged, function norm is NaN\n"); 318184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 319184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 3204b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 321184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 322184914b5SBarry Smith } else if (neP->itflag) { 3234b27c08aSLois Curfman McInnes ierr = SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3241acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 3254b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESConverged_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs); 326184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 327184914b5SBarry Smith } else { 328184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 32952392280SLois Curfman McInnes } 3303a40ed3dSBarry Smith PetscFunctionReturn(0); 33152392280SLois Curfman McInnes } 33252392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 333ebe3b25bSBarry Smith /*MC 334ebe3b25bSBarry Smith SNESTR - Newton based nonlinear solver that uses a trust region 335ebe3b25bSBarry Smith 336ebe3b25bSBarry Smith Options Database: 337ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 338ebe3b25bSBarry Smith . -snes_tr_mu <mu> 339ebe3b25bSBarry Smith . -snes_tr_eta <eta> 340ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 341ebe3b25bSBarry Smith . -snes_tr_delta0 <delta0> 342ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 343ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 344ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 345ebe3b25bSBarry Smith 346ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 347ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 348ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 349ebe3b25bSBarry Smith 350ebe3b25bSBarry Smith This is intended as a model implementation, since it does not 351ebe3b25bSBarry Smith necessarily have many of the bells and whistles of other 352ebe3b25bSBarry Smith implementations. 353ebe3b25bSBarry Smith 354*ee3001cbSBarry Smith Level: intermediate 355*ee3001cbSBarry Smith 356ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance() 357ebe3b25bSBarry Smith 358ebe3b25bSBarry Smith M*/ 359fb2e594dSBarry Smith EXTERN_C_BEGIN 3604a2ae208SSatish Balay #undef __FUNCT__ 3614b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR" 3624b27c08aSLois Curfman McInnes int SNESCreate_TR(SNES snes) 3634800dd8cSBarry Smith { 3644b27c08aSLois Curfman McInnes SNES_TR *neP; 36582502324SSatish Balay int ierr; 3664800dd8cSBarry Smith 3673a40ed3dSBarry Smith PetscFunctionBegin; 3684b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_TR; 3694b27c08aSLois Curfman McInnes snes->solve = SNESSolve_TR; 3704b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_TR; 3714b27c08aSLois Curfman McInnes snes->converged = SNESConverged_TR; 3724b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_TR; 3734b27c08aSLois Curfman McInnes snes->view = SNESView_TR; 3745baf8537SBarry Smith snes->nwork = 0; 375fbe28522SBarry Smith 3764b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_TR,&neP);CHKERRQ(ierr); 3774b27c08aSLois Curfman McInnes PetscLogObjectMemory(snes,sizeof(SNES_TR)); 378fbe28522SBarry Smith snes->data = (void*)neP; 379fbe28522SBarry Smith neP->mu = 0.25; 380fbe28522SBarry Smith neP->eta = 0.75; 381fbe28522SBarry Smith neP->delta = 0.0; 382fbe28522SBarry Smith neP->delta0 = 0.2; 383fbe28522SBarry Smith neP->delta1 = 0.3; 384fbe28522SBarry Smith neP->delta2 = 0.75; 385fbe28522SBarry Smith neP->delta3 = 2.0; 386fbe28522SBarry Smith neP->sigma = 0.0001; 387fbe28522SBarry Smith neP->itflag = 0; 38852392280SLois Curfman McInnes neP->rnorm0 = 0; 38952392280SLois Curfman McInnes neP->ttol = 0; 3903a40ed3dSBarry Smith PetscFunctionReturn(0); 3914800dd8cSBarry Smith } 392fb2e594dSBarry Smith EXTERN_C_END 39382bf6240SBarry Smith 394