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) { 22*71f87433Sdalcinl 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 34df60cc22SBarry Smith /* 354b27c08aSLois Curfman McInnes SNESSolve_TR - Implements Newton's Method with a very simple trust 36ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 374800dd8cSBarry Smith 38ddbbdb52SLois Curfman McInnes 394800dd8cSBarry Smith */ 404a2ae208SSatish Balay #undef __FUNCT__ 414b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR" 426849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes) 434800dd8cSBarry Smith { 444b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR*)snes->data; 456b5873e3SBarry Smith Vec X,F,Y,G,TMP,Ytmp; 466849ba73SBarry Smith PetscErrorCode ierr; 47a7cc72afSBarry Smith PetscInt maxits,i,lits; 48112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 49064f8208SBarry Smith PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1; 5079f36460SBarry Smith PetscScalar cnorm; 51df60cc22SBarry Smith KSP ksp; 52184914b5SBarry Smith SNESConvergedReason reason; 53a7cc72afSBarry Smith PetscTruth conv,breakout = PETSC_FALSE; 544800dd8cSBarry Smith 553a40ed3dSBarry Smith PetscFunctionBegin; 56fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 57fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 5839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 59fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 60fbe28522SBarry Smith G = snes->work[1]; 616b5873e3SBarry Smith Ytmp = snes->work[2]; 624800dd8cSBarry Smith 634c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 644c49b128SBarry Smith snes->iter = 0; 654c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 667e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 674800dd8cSBarry Smith 685334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 69cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 700f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 71fbe28522SBarry Smith snes->norm = fnorm; 720f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 734800dd8cSBarry Smith delta = neP->delta0*fnorm; 744800dd8cSBarry Smith neP->delta = delta; 750462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 7694a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 7794b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 784800dd8cSBarry Smith 7970441072SBarry Smith if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 80b37302c6SLois Curfman McInnes 81d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 82d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 83d034289bSBarry Smith 84a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 857c4af3dcSLois Curfman McInnes ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 867c4af3dcSLois Curfman McInnes if (!conv) { 874b27c08aSLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr); 88ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 897c4af3dcSLois Curfman McInnes } 90df60cc22SBarry Smith 914800dd8cSBarry Smith for (i=0; i<maxits; i++) { 9299a96b7cSMatthew Knepley 9399a96b7cSMatthew Knepley /* Call general purpose update function */ 94e7788613SBarry Smith if (snes->ops->update) { 95e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 9699a96b7cSMatthew Knepley } 9799a96b7cSMatthew Knepley 985334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 9994b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1002deea200SLois Curfman McInnes 1012deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 102*71f87433Sdalcinl ierr = SNES_KSPSolve(snes,ksp,F,Ytmp);CHKERRQ(ierr); 103c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 104329f5518SBarry Smith snes->linear_its += lits; 105ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 106064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 107064f8208SBarry Smith norm1 = nrm; 1086b5873e3SBarry Smith while(1) { 109393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 110064f8208SBarry Smith nrm = norm1; 1116b5873e3SBarry Smith 1122deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 113064f8208SBarry Smith if (nrm >= delta) { 114064f8208SBarry Smith nrm = delta/nrm; 115064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 116064f8208SBarry Smith cnorm = nrm; 117ae15b995SBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr); 1182dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 119064f8208SBarry Smith nrm = gpnorm; 120fbe28522SBarry Smith ynorm = delta; 121fbe28522SBarry Smith } else { 122fbe28522SBarry Smith gpnorm = 0.0; 123ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 124064f8208SBarry Smith ynorm = nrm; 125fbe28522SBarry Smith } 12679f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 12781b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1285334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 129cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1304800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 131fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1324800dd8cSBarry Smith 1334800dd8cSBarry Smith /* Update size of trust region */ 1344800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1354800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1364800dd8cSBarry Smith else delta *= neP->delta3; 137ae15b995SBarry Smith ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr); 138ae15b995SBarry Smith ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr); 1394800dd8cSBarry Smith neP->delta = delta; 1404800dd8cSBarry Smith if (rho > neP->sigma) break; 141ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 1426b5873e3SBarry Smith /* check to see if progress is hopeless */ 143ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 144e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 145184914b5SBarry Smith if (reason) { 14652392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 1475ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 148a7cc72afSBarry Smith breakout = PETSC_TRUE; 149454a90a3SBarry Smith break; 15052392280SLois Curfman McInnes } 15150ffb88aSMatthew Knepley snes->numFailures++; 1524800dd8cSBarry Smith } 1531acabf8cSLois Curfman McInnes if (!breakout) { 1544800dd8cSBarry Smith fnorm = gnorm; 1550f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 156c509a162SBarry Smith snes->iter = i+1; 157fbe28522SBarry Smith snes->norm = fnorm; 1580f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15939e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16039e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 161329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1620462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16394a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1644800dd8cSBarry Smith 1654800dd8cSBarry Smith /* Test for convergence */ 166ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 167e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 168184914b5SBarry Smith if (reason) { 16938442cffSBarry Smith break; 17038442cffSBarry Smith } 1711acabf8cSLois Curfman McInnes } else { 1721acabf8cSLois Curfman McInnes break; 1731acabf8cSLois Curfman McInnes } 17438442cffSBarry Smith } 17538442cffSBarry Smith /* Verify solution is in corect location */ 176e15f7bb5SBarry Smith if (X != snes->vec_sol) { 177393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 178e15f7bb5SBarry Smith } 179e15f7bb5SBarry Smith if (F != snes->vec_func) { 180e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 181e15f7bb5SBarry Smith } 18239e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18339e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18452392280SLois Curfman McInnes if (i == maxits) { 185ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 186184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 18752392280SLois Curfman McInnes } 1880f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 189184914b5SBarry Smith snes->reason = reason; 1900f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1913a40ed3dSBarry Smith PetscFunctionReturn(0); 1924800dd8cSBarry Smith } 1934800dd8cSBarry Smith /*------------------------------------------------------------*/ 1944a2ae208SSatish Balay #undef __FUNCT__ 1954b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR" 1966849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes) 1974800dd8cSBarry Smith { 198dfbe8321SBarry Smith PetscErrorCode ierr; 1993a40ed3dSBarry Smith 2003a40ed3dSBarry Smith PetscFunctionBegin; 201410397dcSLisandro Dalcin if (!snes->work) { 20281b6cf68SLois Curfman McInnes snes->nwork = 4; 203d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 204efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 205410397dcSLisandro Dalcin } 20681b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2073a40ed3dSBarry Smith PetscFunctionReturn(0); 2084800dd8cSBarry Smith } 2094800dd8cSBarry Smith /*------------------------------------------------------------*/ 2104a2ae208SSatish Balay #undef __FUNCT__ 2114b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR" 2126849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes) 2134800dd8cSBarry Smith { 214dfbe8321SBarry Smith PetscErrorCode ierr; 2155baf8537SBarry Smith 2163a40ed3dSBarry Smith PetscFunctionBegin; 2175baf8537SBarry Smith if (snes->nwork) { 2184b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2195baf8537SBarry Smith } 220606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2213a40ed3dSBarry Smith PetscFunctionReturn(0); 2224800dd8cSBarry Smith } 2234800dd8cSBarry Smith /*------------------------------------------------------------*/ 2244800dd8cSBarry Smith 2254a2ae208SSatish Balay #undef __FUNCT__ 2264b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR" 2276849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes) 2284800dd8cSBarry Smith { 2294b27c08aSLois Curfman McInnes SNES_TR *ctx = (SNES_TR *)snes->data; 230dfbe8321SBarry Smith PetscErrorCode ierr; 2314800dd8cSBarry Smith 2323a40ed3dSBarry Smith PetscFunctionBegin; 233b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 23487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 2354b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 2364b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 2374b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 2384b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 2394b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 2404b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 2414b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 242b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2433a40ed3dSBarry Smith PetscFunctionReturn(0); 2444800dd8cSBarry Smith } 2454800dd8cSBarry Smith 2464a2ae208SSatish Balay #undef __FUNCT__ 2474b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR" 2486849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer) 249a935fc98SLois Curfman McInnes { 2504b27c08aSLois Curfman McInnes SNES_TR *tr = (SNES_TR *)snes->data; 251dfbe8321SBarry Smith PetscErrorCode ierr; 25232077d6dSBarry Smith PetscTruth iascii; 253a935fc98SLois Curfman McInnes 2543a40ed3dSBarry Smith PetscFunctionBegin; 25532077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 25632077d6dSBarry Smith if (iascii) { 257a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 258a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2595cd90555SBarry Smith } else { 26079a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 26119bcc07fSBarry Smith } 2623a40ed3dSBarry Smith PetscFunctionReturn(0); 263a935fc98SLois Curfman McInnes } 264a935fc98SLois Curfman McInnes 26552392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2664a2ae208SSatish Balay #undef __FUNCT__ 2674b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR" 2688ac28fe4SSatish Balay /*@C 2694b27c08aSLois Curfman McInnes SNESConverged_TR - Monitors the convergence of the trust region 2704b27c08aSLois Curfman McInnes method SNESTR for solving systems of nonlinear equations (default). 27152392280SLois Curfman McInnes 272c7afd0dbSLois Curfman McInnes Collective on SNES 273c7afd0dbSLois Curfman McInnes 27452392280SLois Curfman McInnes Input Parameters: 275c7afd0dbSLois Curfman McInnes + snes - the SNES context 27652392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 27752392280SLois Curfman McInnes . pnorm - 2-norm of current step 27852392280SLois Curfman McInnes . fnorm - 2-norm of function 279c7afd0dbSLois Curfman McInnes - dummy - unused context 280fee21e36SBarry Smith 281184914b5SBarry Smith Output Parameter: 282184914b5SBarry Smith . reason - one of 28370441072SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 284184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 285184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 286184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 287184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 288184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 289184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 29052392280SLois Curfman McInnes 29152392280SLois Curfman McInnes where 292c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 293c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 294c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 295c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 296c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 297c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 29870441072SBarry Smith . abstol - absolute function norm tolerance, 299c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 300c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 301c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 30252392280SLois Curfman McInnes 30315091d37SBarry Smith Level: intermediate 30415091d37SBarry Smith 30552392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 30652392280SLois Curfman McInnes 30752392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 30852392280SLois Curfman McInnes @*/ 30906ee9f85SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 31052392280SLois Curfman McInnes { 3114b27c08aSLois Curfman McInnes SNES_TR *neP = (SNES_TR *)snes->data; 312dfbe8321SBarry Smith PetscErrorCode ierr; 313ddd12767SBarry Smith 3143a40ed3dSBarry Smith PetscFunctionBegin; 315d252947aSBarry Smith if (fnorm != fnorm) { 316ae15b995SBarry Smith ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 317184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 318184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 319ae15b995SBarry Smith ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr); 320184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 321184914b5SBarry Smith } else if (neP->itflag) { 32206ee9f85SBarry Smith ierr = SNESConverged_LS(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 32343e71028SBarry Smith } else if (snes->nfuncs >= snes->max_funcs) { 324ae15b995SBarry Smith ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 325184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 326184914b5SBarry Smith } else { 327184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 32852392280SLois Curfman McInnes } 3293a40ed3dSBarry Smith PetscFunctionReturn(0); 33052392280SLois Curfman McInnes } 33152392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 332ebe3b25bSBarry Smith /*MC 333ebe3b25bSBarry Smith SNESTR - Newton based nonlinear solver that uses a trust region 334ebe3b25bSBarry Smith 335ebe3b25bSBarry Smith Options Database: 336ebe3b25bSBarry Smith + -snes_trtol <tol> Trust region tolerance 337ebe3b25bSBarry Smith . -snes_tr_mu <mu> 338ebe3b25bSBarry Smith . -snes_tr_eta <eta> 339ebe3b25bSBarry Smith . -snes_tr_sigma <sigma> 340ebe3b25bSBarry Smith . -snes_tr_delta0 <delta0> 341ebe3b25bSBarry Smith . -snes_tr_delta1 <delta1> 342ebe3b25bSBarry Smith . -snes_tr_delta2 <delta2> 343ebe3b25bSBarry Smith - -snes_tr_delta3 <delta3> 344ebe3b25bSBarry Smith 345ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 346ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 347ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 348ebe3b25bSBarry Smith 349ebe3b25bSBarry Smith This is intended as a model implementation, since it does not 350ebe3b25bSBarry Smith necessarily have many of the bells and whistles of other 351ebe3b25bSBarry Smith implementations. 352ebe3b25bSBarry Smith 353ee3001cbSBarry Smith Level: intermediate 354ee3001cbSBarry Smith 355ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance() 356ebe3b25bSBarry Smith 357ebe3b25bSBarry Smith M*/ 358fb2e594dSBarry Smith EXTERN_C_BEGIN 3594a2ae208SSatish Balay #undef __FUNCT__ 3604b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR" 36163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_TR(SNES snes) 3624800dd8cSBarry Smith { 3634b27c08aSLois Curfman McInnes SNES_TR *neP; 364dfbe8321SBarry Smith PetscErrorCode ierr; 3654800dd8cSBarry Smith 3663a40ed3dSBarry Smith PetscFunctionBegin; 367e7788613SBarry Smith snes->ops->setup = SNESSetUp_TR; 368e7788613SBarry Smith snes->ops->solve = SNESSolve_TR; 369e7788613SBarry Smith snes->ops->destroy = SNESDestroy_TR; 370e7788613SBarry Smith snes->ops->converged = SNESConverged_TR; 371e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_TR; 372e7788613SBarry Smith snes->ops->view = SNESView_TR; 3735baf8537SBarry Smith snes->nwork = 0; 374fbe28522SBarry Smith 3754b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_TR,&neP);CHKERRQ(ierr); 37652e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_TR));CHKERRQ(ierr); 377fbe28522SBarry Smith snes->data = (void*)neP; 378fbe28522SBarry Smith neP->mu = 0.25; 379fbe28522SBarry Smith neP->eta = 0.75; 380fbe28522SBarry Smith neP->delta = 0.0; 381fbe28522SBarry Smith neP->delta0 = 0.2; 382fbe28522SBarry Smith neP->delta1 = 0.3; 383fbe28522SBarry Smith neP->delta2 = 0.75; 384fbe28522SBarry Smith neP->delta3 = 2.0; 385fbe28522SBarry Smith neP->sigma = 0.0001; 386ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 38752392280SLois Curfman McInnes neP->rnorm0 = 0; 38852392280SLois Curfman McInnes neP->ttol = 0; 3893a40ed3dSBarry Smith PetscFunctionReturn(0); 3904800dd8cSBarry Smith } 391fb2e594dSBarry Smith EXTERN_C_END 39282bf6240SBarry Smith 393