1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*6831982aSBarry Smith static char vcid[] = "$Id: tr.c,v 1.103 1999/10/13 20:38:30 bsmith Exp bsmith $"; 34800dd8cSBarry Smith #endif 44800dd8cSBarry Smith 570f55243SBarry Smith #include "src/snes/impls/tr/tr.h" /*I "snes.h" I*/ 64800dd8cSBarry Smith 7fbe28522SBarry Smith /* 8df60cc22SBarry Smith This convergence test determines if the two norm of the 9df60cc22SBarry Smith solution lies outside the trust region, if so it halts. 10df60cc22SBarry Smith */ 115615d1e5SSatish Balay #undef __FUNC__ 12*6831982aSBarry Smith #define __FUNC__ "SNES_EQ_TR_KSPConverged_Private" 13*6831982aSBarry Smith int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n, double rnorm, void *ctx) 14df60cc22SBarry Smith { 15df60cc22SBarry Smith SNES snes = (SNES) ctx; 1698120f82SLois Curfman McInnes SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 17*6831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data; 18df60cc22SBarry Smith Vec x; 1998120f82SLois Curfman McInnes double norm; 2098120f82SLois Curfman McInnes int ierr, convinfo; 21df60cc22SBarry Smith 223a40ed3dSBarry Smith PetscFunctionBegin; 2398120f82SLois Curfman McInnes if (snes->ksp_ewconv) { 24a8c6a408SBarry Smith if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Eisenstat-Walker onvergence context not created"); 25a8c6a408SBarry Smith if (n == 0) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);} 2698120f82SLois Curfman McInnes kctx->lresid_last = rnorm; 27df60cc22SBarry Smith } 2898120f82SLois Curfman McInnes convinfo = KSPDefaultConverged(ksp,n,rnorm,ctx); 299b04593eSLois Curfman McInnes if (convinfo) { 30*6831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 313a40ed3dSBarry Smith PetscFunctionReturn(convinfo); 329b04593eSLois Curfman McInnes } 33df60cc22SBarry Smith 34a935fc98SLois Curfman McInnes /* Determine norm of solution */ 3578b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 36cddf8d76SBarry Smith ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 37df60cc22SBarry Smith if (norm >= neP->delta) { 38*6831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 39*6831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n", 40981c4779SBarry Smith neP->delta,norm); 413a40ed3dSBarry Smith PetscFunctionReturn(1); 42df60cc22SBarry Smith } 433a40ed3dSBarry Smith PetscFunctionReturn(0); 44df60cc22SBarry Smith } 4582bf6240SBarry Smith 46df60cc22SBarry Smith /* 47f63b844aSLois Curfman McInnes SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust 48ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 494800dd8cSBarry Smith 504800dd8cSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 514800dd8cSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 52ddbbdb52SLois Curfman McInnes of Mathematical Software", Wayne Cowell, editor. 53ddbbdb52SLois Curfman McInnes 544800dd8cSBarry Smith This is intended as a model implementation, since it does not 554800dd8cSBarry Smith necessarily have many of the bells and whistles of other 564800dd8cSBarry Smith implementations. 574800dd8cSBarry Smith */ 585615d1e5SSatish Balay #undef __FUNC__ 595615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_TR" 60f63b844aSLois Curfman McInnes static int SNESSolve_EQ_TR(SNES snes,int *its) 614800dd8cSBarry Smith { 62*6831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *) snes->data; 636b5873e3SBarry Smith Vec X, F, Y, G, TMP, Ytmp; 640462333dSBarry Smith int maxits, i, ierr, lits, breakout = 0; 65112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 660462333dSBarry Smith double rho, fnorm, gnorm, gpnorm, xnorm, delta,norm,ynorm,norm1; 675334005bSBarry Smith Scalar mone = -1.0,cnorm; 68df60cc22SBarry Smith KSP ksp; 69df60cc22SBarry Smith SLES sles; 70184914b5SBarry Smith SNESConvergedReason reason; 714800dd8cSBarry Smith 723a40ed3dSBarry Smith PetscFunctionBegin; 73fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 74fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 7539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 76fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 77fbe28522SBarry Smith G = snes->work[1]; 786b5873e3SBarry Smith Ytmp = snes->work[2]; 794800dd8cSBarry Smith 807e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 814800dd8cSBarry Smith 825334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 83cddf8d76SBarry Smith ierr = VecNorm(F, NORM_2,&fnorm );CHKERRQ(ierr); /* fnorm <- || F || */ 840f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 85fbe28522SBarry Smith snes->norm = fnorm; 86c509a162SBarry Smith snes->iter = 0; 870f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 884800dd8cSBarry Smith delta = neP->delta0*fnorm; 894800dd8cSBarry Smith neP->delta = delta; 900462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 9194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 924800dd8cSBarry Smith 93184914b5SBarry Smith if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 94b37302c6SLois Curfman McInnes 95d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 96d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 97d034289bSBarry Smith 98a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 9978b31e54SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 10078b31e54SBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 101*6831982aSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 102df60cc22SBarry Smith 1034800dd8cSBarry Smith for ( i=0; i<maxits; i++ ) { 1045334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1055334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1062deea200SLois Curfman McInnes 1072deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10878b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr); 1097a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(lits); 110af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 111cddf8d76SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr); 1122deea200SLois Curfman McInnes norm1 = norm; 1136b5873e3SBarry Smith while(1) { 114393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 1152deea200SLois Curfman McInnes norm = norm1; 1166b5873e3SBarry Smith 1172deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 118eafb4bcbSBarry Smith if (norm >= delta) { 119fbe28522SBarry Smith norm = delta/norm; 120fbe28522SBarry Smith gpnorm = (1.0 - norm)*fnorm; 121eafb4bcbSBarry Smith cnorm = norm; 122af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm ); 123393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 124eafb4bcbSBarry Smith norm = gpnorm; 125fbe28522SBarry Smith ynorm = delta; 126fbe28522SBarry Smith } else { 127fbe28522SBarry Smith gpnorm = 0.0; 128af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n" ); 129fbe28522SBarry Smith ynorm = norm; 130fbe28522SBarry Smith } 1312deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr); /* Y <- X - Y */ 13281b6cf68SLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 1335334005bSBarry Smith ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 134cddf8d76SBarry Smith ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 1354800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 136fbe28522SBarry Smith else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 1374800dd8cSBarry Smith 1384800dd8cSBarry Smith /* Update size of trust region */ 1394800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 1404800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 1414800dd8cSBarry Smith else delta *= neP->delta3; 142af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 143af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1444800dd8cSBarry Smith neP->delta = delta; 1454800dd8cSBarry Smith if (rho > neP->sigma) break; 146af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n"); 1476b5873e3SBarry Smith /* check to see if progress is hopeless */ 14852392280SLois Curfman McInnes neP->itflag = 0; 149184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 150184914b5SBarry Smith if (reason) { 15152392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 152454a90a3SBarry Smith breakout = 1; 153454a90a3SBarry Smith break; 15452392280SLois Curfman McInnes } 15552392280SLois Curfman McInnes snes->nfailures++; 1564800dd8cSBarry Smith } 1571acabf8cSLois Curfman McInnes if (!breakout) { 1584800dd8cSBarry Smith fnorm = gnorm; 1590f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 160c509a162SBarry Smith snes->iter = i+1; 161fbe28522SBarry Smith snes->norm = fnorm; 1620f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16339e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16439e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 165cddf8d76SBarry Smith VecNorm(X, NORM_2,&xnorm ); /* xnorm = || X || */ 1660462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16794a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1684800dd8cSBarry Smith 1694800dd8cSBarry Smith /* Test for convergence */ 17052392280SLois Curfman McInnes neP->itflag = 1; 171184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 172184914b5SBarry Smith if (reason) { 17338442cffSBarry Smith break; 17438442cffSBarry Smith } 1751acabf8cSLois Curfman McInnes } else { 1761acabf8cSLois Curfman McInnes break; 1771acabf8cSLois Curfman McInnes } 17838442cffSBarry Smith } 17939e2f89bSBarry Smith if (X != snes->vec_sol) { 18038442cffSBarry Smith /* Verify solution is in corect location */ 181393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 18239e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18339e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18439e2f89bSBarry Smith } 18552392280SLois Curfman McInnes if (i == maxits) { 186af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits); 18752392280SLois Curfman McInnes i--; 188184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 18952392280SLois Curfman McInnes } 19052392280SLois Curfman McInnes *its = i+1; 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 /*------------------------------------------------------------*/ 1975615d1e5SSatish Balay #undef __FUNC__ 1985615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_TR" 199f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_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); 206393d2d9aSLois Curfman McInnes PLogObjectParents(snes,snes->nwork,snes->work); 20781b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 2094800dd8cSBarry Smith } 2104800dd8cSBarry Smith /*------------------------------------------------------------*/ 2115615d1e5SSatish Balay #undef __FUNC__ 212d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_TR" 213e1311b90SBarry Smith static int SNESDestroy_EQ_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 2265615d1e5SSatish Balay #undef __FUNC__ 2275615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_TR" 228f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes) 2294800dd8cSBarry Smith { 230*6831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 231fbe28522SBarry Smith double tmp; 23225cdf11fSBarry Smith int ierr,flg; 2334800dd8cSBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 23509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_mu",&tmp, &flg);CHKERRQ(ierr); 23625cdf11fSBarry Smith if (flg) {ctx->mu = tmp;} 23709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_eta",&tmp, &flg);CHKERRQ(ierr); 23825cdf11fSBarry Smith if (flg) {ctx->eta = tmp;} 23909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_sigma",&tmp, &flg);CHKERRQ(ierr); 24025cdf11fSBarry Smith if (flg) {ctx->sigma = tmp;} 24109d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta0",&tmp, &flg);CHKERRQ(ierr); 24225cdf11fSBarry Smith if (flg) {ctx->delta0 = tmp;} 24309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta1",&tmp, &flg);CHKERRQ(ierr); 24425cdf11fSBarry Smith if (flg) {ctx->delta1 = tmp;} 24509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta2",&tmp, &flg);CHKERRQ(ierr); 24625cdf11fSBarry Smith if (flg) {ctx->delta2 = tmp;} 24709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta3",&tmp, &flg);CHKERRQ(ierr); 24825cdf11fSBarry Smith if (flg) {ctx->delta3 = tmp;} 2493a40ed3dSBarry Smith PetscFunctionReturn(0); 2504800dd8cSBarry Smith } 2514800dd8cSBarry Smith 2525615d1e5SSatish Balay #undef __FUNC__ 253d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_TR" 254f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_TR(SNES snes,char *p) 2554800dd8cSBarry Smith { 256*6831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 257d132466eSBarry Smith int ierr; 258d132466eSBarry Smith MPI_Comm comm = snes->comm; 2596daaf66cSBarry Smith 2603a40ed3dSBarry Smith PetscFunctionBegin; 261*6831982aSBarry Smith ierr = (*PetscHelpPrintf)(comm," method SNESEQTR (tr) for systems of nonlinear equations:\n");CHKERRQ(ierr); 262d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_mu <mu> (default %g)\n",p,ctx->mu);CHKERRQ(ierr); 263d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_eta <eta> (default %g)\n",p,ctx->eta);CHKERRQ(ierr); 264d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_sigma <sigma> (default %g)\n",p,ctx->sigma);CHKERRQ(ierr); 265d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta0 <delta0> (default %g)\n",p,ctx->delta0);CHKERRQ(ierr); 266d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta1 <delta1> (default %g)\n",p,ctx->delta1);CHKERRQ(ierr); 267d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta2 <delta2> (default %g)\n",p,ctx->delta2);CHKERRQ(ierr); 268d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta3 <delta3> (default %g)\n",p,ctx->delta3);CHKERRQ(ierr); 2693a40ed3dSBarry Smith PetscFunctionReturn(0); 2704800dd8cSBarry Smith } 2714800dd8cSBarry Smith 2725615d1e5SSatish Balay #undef __FUNC__ 273d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_TR" 274e1311b90SBarry Smith static int SNESView_EQ_TR(SNES snes,Viewer viewer) 275a935fc98SLois Curfman McInnes { 276*6831982aSBarry Smith SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data; 27751695f54SLois Curfman McInnes int ierr; 278*6831982aSBarry Smith PetscTruth isascii; 279a935fc98SLois Curfman McInnes 2803a40ed3dSBarry Smith PetscFunctionBegin; 281*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 2820f5bd95cSBarry Smith if (isascii) { 283d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 284d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2855cd90555SBarry Smith } else { 2860f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 28719bcc07fSBarry Smith } 2883a40ed3dSBarry Smith PetscFunctionReturn(0); 289a935fc98SLois Curfman McInnes } 290a935fc98SLois Curfman McInnes 29152392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2925615d1e5SSatish Balay #undef __FUNC__ 2935615d1e5SSatish Balay #define __FUNC__ "SNESConverged_EQ_TR" 29452392280SLois Curfman McInnes /*@ 295f525115eSLois Curfman McInnes SNESConverged_EQ_TR - Monitors the convergence of the trust region 296*6831982aSBarry Smith method SNESEQTR for solving systems of nonlinear equations (default). 29752392280SLois Curfman McInnes 298c7afd0dbSLois Curfman McInnes Collective on SNES 299c7afd0dbSLois Curfman McInnes 30052392280SLois Curfman McInnes Input Parameters: 301c7afd0dbSLois Curfman McInnes + snes - the SNES context 30252392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 30352392280SLois Curfman McInnes . pnorm - 2-norm of current step 30452392280SLois Curfman McInnes . fnorm - 2-norm of function 305c7afd0dbSLois Curfman McInnes - dummy - unused context 306fee21e36SBarry Smith 307184914b5SBarry Smith Output Parameter: 308184914b5SBarry Smith . reason - one of 309184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - ( fnorm < atol ), 310184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - ( pnorm < xtol*xnorm ), 311184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - ( fnorm < rtol*fnorm0 ), 312184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - ( nfct > maxf ), 313184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - ( fnorm == NaN ), 314184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - ( delta < xnorm*deltatol ), 315184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - ( otherwise ) 31652392280SLois Curfman McInnes 31752392280SLois Curfman McInnes where 318c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 319c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 320c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 321c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 322c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 323c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 324c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 325c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 326c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 327c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 32852392280SLois Curfman McInnes 32915091d37SBarry Smith Level: intermediate 33015091d37SBarry Smith 33152392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 33252392280SLois Curfman McInnes 33352392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 33452392280SLois Curfman McInnes @*/ 335184914b5SBarry Smith int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy) 33652392280SLois Curfman McInnes { 337*6831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data; 338184914b5SBarry Smith int ierr; 339ddd12767SBarry Smith 3403a40ed3dSBarry Smith PetscFunctionBegin; 3413a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 342a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3433a40ed3dSBarry Smith } 34452392280SLois Curfman McInnes 345d252947aSBarry Smith if (fnorm != fnorm) { 346981c4779SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n"); 347184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 348184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 349184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 350184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 351184914b5SBarry Smith } else if (neP->itflag) { 352184914b5SBarry Smith ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3531acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 354184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs, snes->max_funcs ); 355184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 356184914b5SBarry Smith } else { 357184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 35852392280SLois Curfman McInnes } 3593a40ed3dSBarry Smith PetscFunctionReturn(0); 36052392280SLois Curfman McInnes } 36152392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 362fb2e594dSBarry Smith EXTERN_C_BEGIN 3635615d1e5SSatish Balay #undef __FUNC__ 3645615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_TR" 365f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes ) 3664800dd8cSBarry Smith { 367*6831982aSBarry Smith SNES_EQ_TR *neP; 3684800dd8cSBarry Smith 3693a40ed3dSBarry Smith PetscFunctionBegin; 3703a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 371a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3723a40ed3dSBarry Smith } 373f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_TR; 374f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_TR; 375f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_TR; 376f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_TR; 377f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_TR; 378f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_TR; 379f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_TR; 3805baf8537SBarry Smith snes->nwork = 0; 381fbe28522SBarry Smith 382*6831982aSBarry Smith neP = PetscNew(SNES_EQ_TR);CHKPTRQ(neP); 383*6831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_TR)); 384fbe28522SBarry Smith snes->data = (void *) neP; 385fbe28522SBarry Smith neP->mu = 0.25; 386fbe28522SBarry Smith neP->eta = 0.75; 387fbe28522SBarry Smith neP->delta = 0.0; 388fbe28522SBarry Smith neP->delta0 = 0.2; 389fbe28522SBarry Smith neP->delta1 = 0.3; 390fbe28522SBarry Smith neP->delta2 = 0.75; 391fbe28522SBarry Smith neP->delta3 = 2.0; 392fbe28522SBarry Smith neP->sigma = 0.0001; 393fbe28522SBarry Smith neP->itflag = 0; 39452392280SLois Curfman McInnes neP->rnorm0 = 0; 39552392280SLois Curfman McInnes neP->ttol = 0; 3963a40ed3dSBarry Smith PetscFunctionReturn(0); 3974800dd8cSBarry Smith } 398fb2e594dSBarry Smith EXTERN_C_END 39982bf6240SBarry Smith 400