1*329f5518SBarry Smith /*$Id: tr.c,v 1.110 1999/12/08 22:17:39 balay Exp bsmith $*/ 24800dd8cSBarry Smith 370f55243SBarry Smith #include "src/snes/impls/tr/tr.h" /*I "snes.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 */ 95615d1e5SSatish Balay #undef __FUNC__ 106831982aSBarry Smith #define __FUNC__ "SNES_EQ_TR_KSPConverged_Private" 113acb151aSSatish Balay int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,double rnorm,KSPConvergedReason *reason,void *ctx) 12df60cc22SBarry Smith { 13df60cc22SBarry Smith SNES snes = (SNES) ctx; 1498120f82SLois Curfman McInnes SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 156831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data; 16df60cc22SBarry Smith Vec x; 1798120f82SLois Curfman McInnes double norm; 183acb151aSSatish Balay int ierr; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 2198120f82SLois Curfman McInnes if (snes->ksp_ewconv) { 22a8c6a408SBarry Smith if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"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 } 26*329f5518SBarry Smith ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr); 27*329f5518SBarry Smith if (*reason) { 28*329f5518SBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%g\n",n,rnorm); 299b04593eSLois Curfman McInnes } 30df60cc22SBarry Smith 31a935fc98SLois Curfman McInnes /* Determine norm of solution */ 3278b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 33cddf8d76SBarry Smith ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 34df60cc22SBarry Smith if (norm >= neP->delta) { 356831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 36c38d4ed2SBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,norm); 37*329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 38df60cc22SBarry Smith } 393a40ed3dSBarry Smith PetscFunctionReturn(0); 40df60cc22SBarry Smith } 4182bf6240SBarry Smith 42df60cc22SBarry Smith /* 43f63b844aSLois Curfman McInnes SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust 44ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 454800dd8cSBarry Smith 464800dd8cSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 474800dd8cSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 48ddbbdb52SLois Curfman McInnes of Mathematical Software", Wayne Cowell, editor. 49ddbbdb52SLois Curfman McInnes 504800dd8cSBarry Smith This is intended as a model implementation, since it does not 514800dd8cSBarry Smith necessarily have many of the bells and whistles of other 524800dd8cSBarry Smith implementations. 534800dd8cSBarry Smith */ 545615d1e5SSatish Balay #undef __FUNC__ 555615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_TR" 56f63b844aSLois Curfman McInnes static int SNESSolve_EQ_TR(SNES snes,int *its) 574800dd8cSBarry Smith { 586831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data; 596b5873e3SBarry Smith Vec X,F,Y,G,TMP,Ytmp; 600462333dSBarry Smith int maxits,i,ierr,lits,breakout = 0; 61112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 620462333dSBarry Smith double rho,fnorm,gnorm,gpnorm,xnorm,delta,norm,ynorm,norm1; 635334005bSBarry Smith Scalar mone = -1.0,cnorm; 64df60cc22SBarry Smith KSP ksp; 65df60cc22SBarry Smith SLES sles; 66184914b5SBarry Smith SNESConvergedReason reason; 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 767e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 774800dd8cSBarry Smith 785334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 79cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 800f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 81fbe28522SBarry Smith snes->norm = fnorm; 82c509a162SBarry Smith snes->iter = 0; 830f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 844800dd8cSBarry Smith delta = neP->delta0*fnorm; 854800dd8cSBarry Smith neP->delta = delta; 860462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 8794a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 884800dd8cSBarry Smith 89184914b5SBarry Smith if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 90b37302c6SLois Curfman McInnes 91d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 92d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 93d034289bSBarry Smith 94a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 9578b31e54SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 9678b31e54SBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 976831982aSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 98df60cc22SBarry Smith 994800dd8cSBarry Smith for (i=0; i<maxits; i++) { 1005334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1015334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1022deea200SLois Curfman McInnes 1032deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10478b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr); 105*329f5518SBarry Smith snes->linear_its += lits; 106af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 107cddf8d76SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr); 1082deea200SLois Curfman McInnes norm1 = norm; 1096b5873e3SBarry Smith while(1) { 110393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 1112deea200SLois Curfman McInnes norm = norm1; 1126b5873e3SBarry Smith 1132deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 114eafb4bcbSBarry Smith if (norm >= delta) { 115fbe28522SBarry Smith norm = delta/norm; 116fbe28522SBarry Smith gpnorm = (1.0 - norm)*fnorm; 117eafb4bcbSBarry Smith cnorm = norm; 118af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm); 119393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 120eafb4bcbSBarry Smith norm = gpnorm; 121fbe28522SBarry Smith ynorm = delta; 122fbe28522SBarry Smith } else { 123fbe28522SBarry Smith gpnorm = 0.0; 124af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n"); 125fbe28522SBarry Smith ynorm = norm; 126fbe28522SBarry Smith } 1272deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr); /* Y <- X - Y */ 12881b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1295334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 130cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1314800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 132fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1334800dd8cSBarry Smith 1344800dd8cSBarry Smith /* Update size of trust region */ 1354800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1364800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1374800dd8cSBarry Smith else delta *= neP->delta3; 138af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 139af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1404800dd8cSBarry Smith neP->delta = delta; 1414800dd8cSBarry Smith if (rho > neP->sigma) break; 142af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n"); 1436b5873e3SBarry Smith /* check to see if progress is hopeless */ 14452392280SLois Curfman McInnes neP->itflag = 0; 145184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 146184914b5SBarry Smith if (reason) { 14752392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 148454a90a3SBarry Smith breakout = 1; 149454a90a3SBarry Smith break; 15052392280SLois Curfman McInnes } 15152392280SLois Curfman McInnes snes->nfailures++; 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; 161*329f5518SBarry 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 */ 16652392280SLois Curfman McInnes neP->itflag = 1; 167184914b5SBarry Smith ierr = (*snes->converged)(snes,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 } 17539e2f89bSBarry Smith if (X != snes->vec_sol) { 17638442cffSBarry Smith /* Verify solution is in corect location */ 177393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 17839e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 17939e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18039e2f89bSBarry Smith } 18152392280SLois Curfman McInnes if (i == maxits) { 182af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits); 18352392280SLois Curfman McInnes i--; 184184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 18552392280SLois Curfman McInnes } 18652392280SLois Curfman McInnes *its = i+1; 1870f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 188184914b5SBarry Smith snes->reason = reason; 1890f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 1914800dd8cSBarry Smith } 1924800dd8cSBarry Smith /*------------------------------------------------------------*/ 1935615d1e5SSatish Balay #undef __FUNC__ 1945615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_TR" 195f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_TR(SNES snes) 1964800dd8cSBarry Smith { 197fbe28522SBarry Smith int ierr; 1983a40ed3dSBarry Smith 1993a40ed3dSBarry Smith PetscFunctionBegin; 20081b6cf68SLois Curfman McInnes snes->nwork = 4; 201d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 202393d2d9aSLois Curfman McInnes PLogObjectParents(snes,snes->nwork,snes->work); 20381b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2043a40ed3dSBarry Smith PetscFunctionReturn(0); 2054800dd8cSBarry Smith } 2064800dd8cSBarry Smith /*------------------------------------------------------------*/ 2075615d1e5SSatish Balay #undef __FUNC__ 208d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_TR" 209e1311b90SBarry Smith static int SNESDestroy_EQ_TR(SNES snes) 2104800dd8cSBarry Smith { 211393d2d9aSLois Curfman McInnes int ierr; 2125baf8537SBarry Smith 2133a40ed3dSBarry Smith PetscFunctionBegin; 2145baf8537SBarry Smith if (snes->nwork) { 2154b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2165baf8537SBarry Smith } 217606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2183a40ed3dSBarry Smith PetscFunctionReturn(0); 2194800dd8cSBarry Smith } 2204800dd8cSBarry Smith /*------------------------------------------------------------*/ 2214800dd8cSBarry Smith 2225615d1e5SSatish Balay #undef __FUNC__ 2235615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_TR" 224f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes) 2254800dd8cSBarry Smith { 2266831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 227fbe28522SBarry Smith double tmp; 228f1af5d2fSBarry Smith int ierr; 229f1af5d2fSBarry Smith PetscTruth flg; 2304800dd8cSBarry Smith 2313a40ed3dSBarry Smith PetscFunctionBegin; 23209d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_mu",&tmp,&flg);CHKERRQ(ierr); 23325cdf11fSBarry Smith if (flg) {ctx->mu = tmp;} 23409d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_eta",&tmp,&flg);CHKERRQ(ierr); 23525cdf11fSBarry Smith if (flg) {ctx->eta = tmp;} 23609d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_sigma",&tmp,&flg);CHKERRQ(ierr); 23725cdf11fSBarry Smith if (flg) {ctx->sigma = tmp;} 23809d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta0",&tmp,&flg);CHKERRQ(ierr); 23925cdf11fSBarry Smith if (flg) {ctx->delta0 = tmp;} 24009d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta1",&tmp,&flg);CHKERRQ(ierr); 24125cdf11fSBarry Smith if (flg) {ctx->delta1 = tmp;} 24209d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta2",&tmp,&flg);CHKERRQ(ierr); 24325cdf11fSBarry Smith if (flg) {ctx->delta2 = tmp;} 24409d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta3",&tmp,&flg);CHKERRQ(ierr); 24525cdf11fSBarry Smith if (flg) {ctx->delta3 = tmp;} 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2474800dd8cSBarry Smith } 2484800dd8cSBarry Smith 2495615d1e5SSatish Balay #undef __FUNC__ 250d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_TR" 251f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_TR(SNES snes,char *p) 2524800dd8cSBarry Smith { 2536831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 254d132466eSBarry Smith int ierr; 255d132466eSBarry Smith MPI_Comm comm = snes->comm; 2566daaf66cSBarry Smith 2573a40ed3dSBarry Smith PetscFunctionBegin; 2586831982aSBarry Smith ierr = (*PetscHelpPrintf)(comm," method SNESEQTR (tr) for systems of nonlinear equations:\n");CHKERRQ(ierr); 259d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_mu <mu> (default %g)\n",p,ctx->mu);CHKERRQ(ierr); 260d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_eta <eta> (default %g)\n",p,ctx->eta);CHKERRQ(ierr); 261d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_sigma <sigma> (default %g)\n",p,ctx->sigma);CHKERRQ(ierr); 262d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta0 <delta0> (default %g)\n",p,ctx->delta0);CHKERRQ(ierr); 263d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta1 <delta1> (default %g)\n",p,ctx->delta1);CHKERRQ(ierr); 264d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta2 <delta2> (default %g)\n",p,ctx->delta2);CHKERRQ(ierr); 265d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta3 <delta3> (default %g)\n",p,ctx->delta3);CHKERRQ(ierr); 2663a40ed3dSBarry Smith PetscFunctionReturn(0); 2674800dd8cSBarry Smith } 2684800dd8cSBarry Smith 2695615d1e5SSatish Balay #undef __FUNC__ 270d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_TR" 271e1311b90SBarry Smith static int SNESView_EQ_TR(SNES snes,Viewer viewer) 272a935fc98SLois Curfman McInnes { 2736831982aSBarry Smith SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data; 27451695f54SLois Curfman McInnes int ierr; 2756831982aSBarry Smith PetscTruth isascii; 276a935fc98SLois Curfman McInnes 2773a40ed3dSBarry Smith PetscFunctionBegin; 2786831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 2790f5bd95cSBarry Smith if (isascii) { 280d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 281d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2825cd90555SBarry Smith } else { 2830f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 28419bcc07fSBarry Smith } 2853a40ed3dSBarry Smith PetscFunctionReturn(0); 286a935fc98SLois Curfman McInnes } 287a935fc98SLois Curfman McInnes 28852392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2895615d1e5SSatish Balay #undef __FUNC__ 2905615d1e5SSatish Balay #define __FUNC__ "SNESConverged_EQ_TR" 29152392280SLois Curfman McInnes /*@ 292f525115eSLois Curfman McInnes SNESConverged_EQ_TR - Monitors the convergence of the trust region 2936831982aSBarry Smith method SNESEQTR for solving systems of nonlinear equations (default). 29452392280SLois Curfman McInnes 295c7afd0dbSLois Curfman McInnes Collective on SNES 296c7afd0dbSLois Curfman McInnes 29752392280SLois Curfman McInnes Input Parameters: 298c7afd0dbSLois Curfman McInnes + snes - the SNES context 29952392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 30052392280SLois Curfman McInnes . pnorm - 2-norm of current step 30152392280SLois Curfman McInnes . fnorm - 2-norm of function 302c7afd0dbSLois Curfman McInnes - dummy - unused context 303fee21e36SBarry Smith 304184914b5SBarry Smith Output Parameter: 305184914b5SBarry Smith . reason - one of 306184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol), 307184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 308184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 309184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 310184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 311184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 312184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - (otherwise) 31352392280SLois Curfman McInnes 31452392280SLois Curfman McInnes where 315c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 316c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 317c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 318c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 319c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 320c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 321c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 322c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 323c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 324c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 32552392280SLois Curfman McInnes 32615091d37SBarry Smith Level: intermediate 32715091d37SBarry Smith 32852392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 32952392280SLois Curfman McInnes 33052392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 33152392280SLois Curfman McInnes @*/ 332184914b5SBarry Smith int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy) 33352392280SLois Curfman McInnes { 3346831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data; 335184914b5SBarry Smith int ierr; 336ddd12767SBarry Smith 3373a40ed3dSBarry Smith PetscFunctionBegin; 3383a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 339a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3403a40ed3dSBarry Smith } 34152392280SLois Curfman McInnes 342d252947aSBarry Smith if (fnorm != fnorm) { 343981c4779SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n"); 344184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 345184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 346184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 347184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 348184914b5SBarry Smith } else if (neP->itflag) { 349184914b5SBarry Smith ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3501acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 351184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs); 352184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 353184914b5SBarry Smith } else { 354184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 35552392280SLois Curfman McInnes } 3563a40ed3dSBarry Smith PetscFunctionReturn(0); 35752392280SLois Curfman McInnes } 35852392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 359fb2e594dSBarry Smith EXTERN_C_BEGIN 3605615d1e5SSatish Balay #undef __FUNC__ 3615615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_TR" 362f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes) 3634800dd8cSBarry Smith { 3646831982aSBarry Smith SNES_EQ_TR *neP; 3654800dd8cSBarry Smith 3663a40ed3dSBarry Smith PetscFunctionBegin; 3673a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 368a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3693a40ed3dSBarry Smith } 370f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_TR; 371f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_TR; 372f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_TR; 373f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_TR; 374f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_TR; 375f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_TR; 376f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_TR; 3775baf8537SBarry Smith snes->nwork = 0; 378fbe28522SBarry Smith 3796831982aSBarry Smith neP = PetscNew(SNES_EQ_TR);CHKPTRQ(neP); 3806831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_TR)); 381fbe28522SBarry Smith snes->data = (void*)neP; 382fbe28522SBarry Smith neP->mu = 0.25; 383fbe28522SBarry Smith neP->eta = 0.75; 384fbe28522SBarry Smith neP->delta = 0.0; 385fbe28522SBarry Smith neP->delta0 = 0.2; 386fbe28522SBarry Smith neP->delta1 = 0.3; 387fbe28522SBarry Smith neP->delta2 = 0.75; 388fbe28522SBarry Smith neP->delta3 = 2.0; 389fbe28522SBarry Smith neP->sigma = 0.0001; 390fbe28522SBarry Smith neP->itflag = 0; 39152392280SLois Curfman McInnes neP->rnorm0 = 0; 39252392280SLois Curfman McInnes neP->ttol = 0; 3933a40ed3dSBarry Smith PetscFunctionReturn(0); 3944800dd8cSBarry Smith } 395fb2e594dSBarry Smith EXTERN_C_END 39682bf6240SBarry Smith 397