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; 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; 18dfbe8321SBarry Smith PetscErrorCode ierr; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 2198120f82SLois Curfman McInnes if (snes->ksp_ewconv) { 2290faec13SBarry Smith if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker convergence 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) { 28*ae15b995SBarry Smith ierr = PetscInfo2(snes,"regular convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr); 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) { 35*ae15b995SBarry Smith ierr = PetscInfo2(snes,"KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr); 36*ae15b995SBarry Smith ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr); 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" 506849ba73SBarry Smith static PetscErrorCode 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; 546849ba73SBarry Smith PetscErrorCode ierr; 55a7cc72afSBarry Smith PetscInt maxits,i,lits; 56112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 57064f8208SBarry Smith PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1; 5879f36460SBarry Smith PetscScalar cnorm; 59df60cc22SBarry Smith KSP ksp; 60184914b5SBarry Smith SNESConvergedReason reason; 61a7cc72afSBarry Smith PetscTruth conv,breakout = PETSC_FALSE; 624800dd8cSBarry Smith 633a40ed3dSBarry Smith PetscFunctionBegin; 64fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 65fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 6639e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 67fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 68fbe28522SBarry Smith G = snes->work[1]; 696b5873e3SBarry Smith Ytmp = snes->work[2]; 704800dd8cSBarry Smith 714c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 724c49b128SBarry Smith snes->iter = 0; 734c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 747e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 754800dd8cSBarry Smith 765334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 77cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 780f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 79fbe28522SBarry Smith snes->norm = fnorm; 800f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 814800dd8cSBarry Smith delta = neP->delta0*fnorm; 824800dd8cSBarry Smith neP->delta = delta; 830462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 8494a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 8594b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 864800dd8cSBarry Smith 8770441072SBarry Smith if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 88b37302c6SLois Curfman McInnes 89d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 90d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 91d034289bSBarry Smith 92a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 937c4af3dcSLois Curfman McInnes ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 947c4af3dcSLois Curfman McInnes if (!conv) { 954b27c08aSLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr); 96*ae15b995SBarry Smith ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 977c4af3dcSLois Curfman McInnes } 98df60cc22SBarry Smith 994800dd8cSBarry Smith for (i=0; i<maxits; i++) { 1005334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 10194b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1022deea200SLois Curfman McInnes 1032deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10423ce1328SBarry Smith ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 105c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 106329f5518SBarry Smith snes->linear_its += lits; 107*ae15b995SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 108064f8208SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 109064f8208SBarry Smith norm1 = nrm; 1106b5873e3SBarry Smith while(1) { 111393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 112064f8208SBarry Smith nrm = norm1; 1136b5873e3SBarry Smith 1142deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 115064f8208SBarry Smith if (nrm >= delta) { 116064f8208SBarry Smith nrm = delta/nrm; 117064f8208SBarry Smith gpnorm = (1.0 - nrm)*fnorm; 118064f8208SBarry Smith cnorm = nrm; 119*ae15b995SBarry Smith ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr); 1202dcb1b2aSMatthew Knepley ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 121064f8208SBarry Smith nrm = gpnorm; 122fbe28522SBarry Smith ynorm = delta; 123fbe28522SBarry Smith } else { 124fbe28522SBarry Smith gpnorm = 0.0; 125*ae15b995SBarry Smith ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 126064f8208SBarry Smith ynorm = nrm; 127fbe28522SBarry Smith } 12879f36460SBarry Smith ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 12981b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1305334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 131cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1324800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 133fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1344800dd8cSBarry Smith 1354800dd8cSBarry Smith /* Update size of trust region */ 1364800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1374800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1384800dd8cSBarry Smith else delta *= neP->delta3; 139*ae15b995SBarry Smith ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr); 140*ae15b995SBarry Smith ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr); 1414800dd8cSBarry Smith neP->delta = delta; 1424800dd8cSBarry Smith if (rho > neP->sigma) break; 143*ae15b995SBarry Smith ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 1446b5873e3SBarry Smith /* check to see if progress is hopeless */ 145ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 146184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 147184914b5SBarry Smith if (reason) { 14852392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 1495ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 150a7cc72afSBarry Smith breakout = PETSC_TRUE; 151454a90a3SBarry Smith break; 15252392280SLois Curfman McInnes } 15350ffb88aSMatthew Knepley snes->numFailures++; 1544800dd8cSBarry Smith } 1551acabf8cSLois Curfman McInnes if (!breakout) { 1564800dd8cSBarry Smith fnorm = gnorm; 1570f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 158c509a162SBarry Smith snes->iter = i+1; 159fbe28522SBarry Smith snes->norm = fnorm; 1600f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16139e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16239e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 163329f5518SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1640462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16594a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1664800dd8cSBarry Smith 1674800dd8cSBarry Smith /* Test for convergence */ 168ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 169184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 170184914b5SBarry Smith if (reason) { 17138442cffSBarry Smith break; 17238442cffSBarry Smith } 1731acabf8cSLois Curfman McInnes } else { 1741acabf8cSLois Curfman McInnes break; 1751acabf8cSLois Curfman McInnes } 17638442cffSBarry Smith } 17738442cffSBarry Smith /* Verify solution is in corect location */ 178e15f7bb5SBarry Smith if (X != snes->vec_sol) { 179393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 180e15f7bb5SBarry Smith } 181e15f7bb5SBarry Smith if (F != snes->vec_func) { 182e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 183e15f7bb5SBarry Smith } 18439e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18539e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18652392280SLois Curfman McInnes if (i == maxits) { 187*ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 188184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 18952392280SLois Curfman McInnes } 1900f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 191184914b5SBarry Smith snes->reason = reason; 1920f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 1944800dd8cSBarry Smith } 1954800dd8cSBarry Smith /*------------------------------------------------------------*/ 1964a2ae208SSatish Balay #undef __FUNCT__ 1974b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR" 1986849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes) 1994800dd8cSBarry Smith { 200dfbe8321SBarry Smith PetscErrorCode ierr; 2013a40ed3dSBarry Smith 2023a40ed3dSBarry Smith PetscFunctionBegin; 20381b6cf68SLois Curfman McInnes snes->nwork = 4; 204d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 205efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 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 @*/ 30963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,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) { 316*ae15b995SBarry 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) { 319*ae15b995SBarry 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) { 3224b27c08aSLois Curfman McInnes ierr = SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 32343e71028SBarry Smith } else if (snes->nfuncs >= snes->max_funcs) { 324*ae15b995SBarry 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; 3674b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_TR; 3684b27c08aSLois Curfman McInnes snes->solve = SNESSolve_TR; 3694b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_TR; 3704b27c08aSLois Curfman McInnes snes->converged = SNESConverged_TR; 3714b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_TR; 3724b27c08aSLois Curfman McInnes snes->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