1*c38d4ed2SBarry Smith /*$Id: tr.c,v 1.106 1999/11/05 14:47:12 bsmith 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" 116831982aSBarry Smith int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n, double rnorm, 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; 1898120f82SLois Curfman McInnes int ierr, convinfo; 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"); 23a8c6a408SBarry Smith if (n == 0) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);} 2498120f82SLois Curfman McInnes kctx->lresid_last = rnorm; 25df60cc22SBarry Smith } 2698120f82SLois Curfman McInnes convinfo = KSPDefaultConverged(ksp,n,rnorm,ctx); 279b04593eSLois Curfman McInnes if (convinfo) { 286831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 293a40ed3dSBarry Smith PetscFunctionReturn(convinfo); 309b04593eSLois Curfman McInnes } 31df60cc22SBarry Smith 32a935fc98SLois Curfman McInnes /* Determine norm of solution */ 3378b31e54SBarry Smith ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 34cddf8d76SBarry Smith ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 35df60cc22SBarry Smith if (norm >= neP->delta) { 366831982aSBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 37*c38d4ed2SBarry Smith PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,norm); 383a40ed3dSBarry Smith PetscFunctionReturn(1); 39df60cc22SBarry Smith } 403a40ed3dSBarry Smith PetscFunctionReturn(0); 41df60cc22SBarry Smith } 4282bf6240SBarry Smith 43df60cc22SBarry Smith /* 44f63b844aSLois Curfman McInnes SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust 45ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 464800dd8cSBarry Smith 474800dd8cSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 484800dd8cSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 49ddbbdb52SLois Curfman McInnes of Mathematical Software", Wayne Cowell, editor. 50ddbbdb52SLois Curfman McInnes 514800dd8cSBarry Smith This is intended as a model implementation, since it does not 524800dd8cSBarry Smith necessarily have many of the bells and whistles of other 534800dd8cSBarry Smith implementations. 544800dd8cSBarry Smith */ 555615d1e5SSatish Balay #undef __FUNC__ 565615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_TR" 57f63b844aSLois Curfman McInnes static int SNESSolve_EQ_TR(SNES snes,int *its) 584800dd8cSBarry Smith { 596831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *) snes->data; 606b5873e3SBarry Smith Vec X, F, Y, G, TMP, Ytmp; 610462333dSBarry Smith int maxits, i, ierr, lits, breakout = 0; 62112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 630462333dSBarry Smith double rho, fnorm, gnorm, gpnorm, xnorm, delta,norm,ynorm,norm1; 645334005bSBarry Smith Scalar mone = -1.0,cnorm; 65df60cc22SBarry Smith KSP ksp; 66df60cc22SBarry Smith SLES sles; 67184914b5SBarry Smith SNESConvergedReason reason; 684800dd8cSBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 70fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 71fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 7239e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 73fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 74fbe28522SBarry Smith G = snes->work[1]; 756b5873e3SBarry Smith Ytmp = snes->work[2]; 764800dd8cSBarry Smith 777e6560a3SBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 784800dd8cSBarry Smith 795334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 80cddf8d76SBarry Smith ierr = VecNorm(F, NORM_2,&fnorm );CHKERRQ(ierr); /* fnorm <- || F || */ 810f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 82fbe28522SBarry Smith snes->norm = fnorm; 83c509a162SBarry Smith snes->iter = 0; 840f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 854800dd8cSBarry Smith delta = neP->delta0*fnorm; 864800dd8cSBarry Smith neP->delta = delta; 870462333dSBarry Smith SNESLogConvHistory(snes,fnorm,0); 8894a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 894800dd8cSBarry Smith 90184914b5SBarry Smith if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 91b37302c6SLois Curfman McInnes 92d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 93d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 94d034289bSBarry Smith 95a935fc98SLois Curfman McInnes /* Set the stopping criteria to use the More' trick. */ 9678b31e54SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 9778b31e54SBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 986831982aSBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr); 99df60cc22SBarry Smith 1004800dd8cSBarry Smith for ( i=0; i<maxits; i++ ) { 1015334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1025334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 1032deea200SLois Curfman McInnes 1042deea200SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 10578b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr); 1067a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(lits); 107af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 108cddf8d76SBarry Smith ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr); 1092deea200SLois Curfman McInnes norm1 = norm; 1106b5873e3SBarry Smith while(1) { 111393d2d9aSLois Curfman McInnes ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 1122deea200SLois Curfman McInnes norm = norm1; 1136b5873e3SBarry Smith 1142deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 115eafb4bcbSBarry Smith if (norm >= delta) { 116fbe28522SBarry Smith norm = delta/norm; 117fbe28522SBarry Smith gpnorm = (1.0 - norm)*fnorm; 118eafb4bcbSBarry Smith cnorm = norm; 119af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm ); 120393d2d9aSLois Curfman McInnes ierr = VecScale(&cnorm,Y);CHKERRQ(ierr); 121eafb4bcbSBarry Smith norm = gpnorm; 122fbe28522SBarry Smith ynorm = delta; 123fbe28522SBarry Smith } else { 124fbe28522SBarry Smith gpnorm = 0.0; 125af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n" ); 126fbe28522SBarry Smith ynorm = norm; 127fbe28522SBarry Smith } 1282deea200SLois Curfman McInnes ierr = VecAYPX(&mone,X,Y);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; 139af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm); 140af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta); 1414800dd8cSBarry Smith neP->delta = delta; 1424800dd8cSBarry Smith if (rho > neP->sigma) break; 143af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n"); 1446b5873e3SBarry Smith /* check to see if progress is hopeless */ 14552392280SLois Curfman McInnes neP->itflag = 0; 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 */ 149454a90a3SBarry Smith breakout = 1; 150454a90a3SBarry Smith break; 15152392280SLois Curfman McInnes } 15252392280SLois Curfman McInnes snes->nfailures++; 1534800dd8cSBarry Smith } 1541acabf8cSLois Curfman McInnes if (!breakout) { 1554800dd8cSBarry Smith fnorm = gnorm; 1560f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 157c509a162SBarry Smith snes->iter = i+1; 158fbe28522SBarry Smith snes->norm = fnorm; 1590f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16039e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16139e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 162cddf8d76SBarry Smith VecNorm(X, NORM_2,&xnorm ); /* xnorm = || X || */ 1630462333dSBarry Smith SNESLogConvHistory(snes,fnorm,lits); 16494a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1654800dd8cSBarry Smith 1664800dd8cSBarry Smith /* Test for convergence */ 16752392280SLois Curfman McInnes neP->itflag = 1; 168184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 169184914b5SBarry Smith if (reason) { 17038442cffSBarry Smith break; 17138442cffSBarry Smith } 1721acabf8cSLois Curfman McInnes } else { 1731acabf8cSLois Curfman McInnes break; 1741acabf8cSLois Curfman McInnes } 17538442cffSBarry Smith } 17639e2f89bSBarry Smith if (X != snes->vec_sol) { 17738442cffSBarry Smith /* Verify solution is in corect location */ 178393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 17939e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 18039e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 18139e2f89bSBarry Smith } 18252392280SLois Curfman McInnes if (i == maxits) { 183af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits); 18452392280SLois Curfman McInnes i--; 185184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 18652392280SLois Curfman McInnes } 18752392280SLois Curfman McInnes *its = i+1; 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 /*------------------------------------------------------------*/ 1945615d1e5SSatish Balay #undef __FUNC__ 1955615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_TR" 196f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_TR(SNES snes) 1974800dd8cSBarry Smith { 198fbe28522SBarry Smith int ierr; 1993a40ed3dSBarry Smith 2003a40ed3dSBarry Smith PetscFunctionBegin; 20181b6cf68SLois Curfman McInnes snes->nwork = 4; 202d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work );CHKERRQ(ierr); 203393d2d9aSLois Curfman McInnes PLogObjectParents(snes,snes->nwork,snes->work); 20481b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2053a40ed3dSBarry Smith PetscFunctionReturn(0); 2064800dd8cSBarry Smith } 2074800dd8cSBarry Smith /*------------------------------------------------------------*/ 2085615d1e5SSatish Balay #undef __FUNC__ 209d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_TR" 210e1311b90SBarry Smith static int SNESDestroy_EQ_TR(SNES snes ) 2114800dd8cSBarry Smith { 212393d2d9aSLois Curfman McInnes int ierr; 2135baf8537SBarry Smith 2143a40ed3dSBarry Smith PetscFunctionBegin; 2155baf8537SBarry Smith if (snes->nwork) { 2164b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 2175baf8537SBarry Smith } 218606d414cSSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2193a40ed3dSBarry Smith PetscFunctionReturn(0); 2204800dd8cSBarry Smith } 2214800dd8cSBarry Smith /*------------------------------------------------------------*/ 2224800dd8cSBarry Smith 2235615d1e5SSatish Balay #undef __FUNC__ 2245615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_TR" 225f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes) 2264800dd8cSBarry Smith { 2276831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 228fbe28522SBarry Smith double tmp; 229f1af5d2fSBarry Smith int ierr; 230f1af5d2fSBarry Smith PetscTruth flg; 2314800dd8cSBarry Smith 2323a40ed3dSBarry Smith PetscFunctionBegin; 23309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_mu",&tmp, &flg);CHKERRQ(ierr); 23425cdf11fSBarry Smith if (flg) {ctx->mu = tmp;} 23509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_eta",&tmp, &flg);CHKERRQ(ierr); 23625cdf11fSBarry Smith if (flg) {ctx->eta = tmp;} 23709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_sigma",&tmp, &flg);CHKERRQ(ierr); 23825cdf11fSBarry Smith if (flg) {ctx->sigma = tmp;} 23909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta0",&tmp, &flg);CHKERRQ(ierr); 24025cdf11fSBarry Smith if (flg) {ctx->delta0 = tmp;} 24109d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta1",&tmp, &flg);CHKERRQ(ierr); 24225cdf11fSBarry Smith if (flg) {ctx->delta1 = tmp;} 24309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta2",&tmp, &flg);CHKERRQ(ierr); 24425cdf11fSBarry Smith if (flg) {ctx->delta2 = tmp;} 24509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_tr_delta3",&tmp, &flg);CHKERRQ(ierr); 24625cdf11fSBarry Smith if (flg) {ctx->delta3 = tmp;} 2473a40ed3dSBarry Smith PetscFunctionReturn(0); 2484800dd8cSBarry Smith } 2494800dd8cSBarry Smith 2505615d1e5SSatish Balay #undef __FUNC__ 251d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_TR" 252f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_TR(SNES snes,char *p) 2534800dd8cSBarry Smith { 2546831982aSBarry Smith SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data; 255d132466eSBarry Smith int ierr; 256d132466eSBarry Smith MPI_Comm comm = snes->comm; 2576daaf66cSBarry Smith 2583a40ed3dSBarry Smith PetscFunctionBegin; 2596831982aSBarry Smith ierr = (*PetscHelpPrintf)(comm," method SNESEQTR (tr) for systems of nonlinear equations:\n");CHKERRQ(ierr); 260d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_mu <mu> (default %g)\n",p,ctx->mu);CHKERRQ(ierr); 261d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_eta <eta> (default %g)\n",p,ctx->eta);CHKERRQ(ierr); 262d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_sigma <sigma> (default %g)\n",p,ctx->sigma);CHKERRQ(ierr); 263d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta0 <delta0> (default %g)\n",p,ctx->delta0);CHKERRQ(ierr); 264d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta1 <delta1> (default %g)\n",p,ctx->delta1);CHKERRQ(ierr); 265d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta2 <delta2> (default %g)\n",p,ctx->delta2);CHKERRQ(ierr); 266d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," %ssnes_eq_tr_delta3 <delta3> (default %g)\n",p,ctx->delta3);CHKERRQ(ierr); 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2684800dd8cSBarry Smith } 2694800dd8cSBarry Smith 2705615d1e5SSatish Balay #undef __FUNC__ 271d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_TR" 272e1311b90SBarry Smith static int SNESView_EQ_TR(SNES snes,Viewer viewer) 273a935fc98SLois Curfman McInnes { 2746831982aSBarry Smith SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data; 27551695f54SLois Curfman McInnes int ierr; 2766831982aSBarry Smith PetscTruth isascii; 277a935fc98SLois Curfman McInnes 2783a40ed3dSBarry Smith PetscFunctionBegin; 2796831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 2800f5bd95cSBarry Smith if (isascii) { 281d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 282d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 2835cd90555SBarry Smith } else { 2840f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 28519bcc07fSBarry Smith } 2863a40ed3dSBarry Smith PetscFunctionReturn(0); 287a935fc98SLois Curfman McInnes } 288a935fc98SLois Curfman McInnes 28952392280SLois Curfman McInnes /* ---------------------------------------------------------------- */ 2905615d1e5SSatish Balay #undef __FUNC__ 2915615d1e5SSatish Balay #define __FUNC__ "SNESConverged_EQ_TR" 29252392280SLois Curfman McInnes /*@ 293f525115eSLois Curfman McInnes SNESConverged_EQ_TR - Monitors the convergence of the trust region 2946831982aSBarry Smith method SNESEQTR for solving systems of nonlinear equations (default). 29552392280SLois Curfman McInnes 296c7afd0dbSLois Curfman McInnes Collective on SNES 297c7afd0dbSLois Curfman McInnes 29852392280SLois Curfman McInnes Input Parameters: 299c7afd0dbSLois Curfman McInnes + snes - the SNES context 30052392280SLois Curfman McInnes . xnorm - 2-norm of current iterate 30152392280SLois Curfman McInnes . pnorm - 2-norm of current step 30252392280SLois Curfman McInnes . fnorm - 2-norm of function 303c7afd0dbSLois Curfman McInnes - dummy - unused context 304fee21e36SBarry Smith 305184914b5SBarry Smith Output Parameter: 306184914b5SBarry Smith . reason - one of 307184914b5SBarry Smith $ SNES_CONVERGED_FNORM_ABS - ( fnorm < atol ), 308184914b5SBarry Smith $ SNES_CONVERGED_PNORM_RELATIVE - ( pnorm < xtol*xnorm ), 309184914b5SBarry Smith $ SNES_CONVERGED_FNORM_RELATIVE - ( fnorm < rtol*fnorm0 ), 310184914b5SBarry Smith $ SNES_DIVERGED_FUNCTION_COUNT - ( nfct > maxf ), 311184914b5SBarry Smith $ SNES_DIVERGED_FNORM_NAN - ( fnorm == NaN ), 312184914b5SBarry Smith $ SNES_CONVERGED_TR_DELTA - ( delta < xnorm*deltatol ), 313184914b5SBarry Smith $ SNES_CONVERGED_ITERATING - ( otherwise ) 31452392280SLois Curfman McInnes 31552392280SLois Curfman McInnes where 316c7afd0dbSLois Curfman McInnes + delta - trust region paramenter 317c7afd0dbSLois Curfman McInnes . deltatol - trust region size tolerance, 318c7afd0dbSLois Curfman McInnes set with SNESSetTrustRegionTolerance() 319c7afd0dbSLois Curfman McInnes . maxf - maximum number of function evaluations, 320c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 321c7afd0dbSLois Curfman McInnes . nfct - number of function evaluations, 322c7afd0dbSLois Curfman McInnes . atol - absolute function norm tolerance, 323c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 324c7afd0dbSLois Curfman McInnes - xtol - relative function norm tolerance, 325c7afd0dbSLois Curfman McInnes set with SNESSetTolerances() 32652392280SLois Curfman McInnes 32715091d37SBarry Smith Level: intermediate 32815091d37SBarry Smith 32952392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence 33052392280SLois Curfman McInnes 33152392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 33252392280SLois Curfman McInnes @*/ 333184914b5SBarry Smith int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy) 33452392280SLois Curfman McInnes { 3356831982aSBarry Smith SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data; 336184914b5SBarry Smith int ierr; 337ddd12767SBarry Smith 3383a40ed3dSBarry Smith PetscFunctionBegin; 3393a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 340a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3413a40ed3dSBarry Smith } 34252392280SLois Curfman McInnes 343d252947aSBarry Smith if (fnorm != fnorm) { 344981c4779SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n"); 345184914b5SBarry Smith *reason = SNES_DIVERGED_FNORM_NAN; 346184914b5SBarry Smith } else if (neP->delta < xnorm * snes->deltatol) { 347184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol); 348184914b5SBarry Smith *reason = SNES_CONVERGED_TR_DELTA; 349184914b5SBarry Smith } else if (neP->itflag) { 350184914b5SBarry Smith ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 3511acabf8cSLois Curfman McInnes } else if (snes->nfuncs > snes->max_funcs) { 352184914b5SBarry Smith PLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs, snes->max_funcs ); 353184914b5SBarry Smith *reason = SNES_DIVERGED_FUNCTION_COUNT; 354184914b5SBarry Smith } else { 355184914b5SBarry Smith *reason = SNES_CONVERGED_ITERATING; 35652392280SLois Curfman McInnes } 3573a40ed3dSBarry Smith PetscFunctionReturn(0); 35852392280SLois Curfman McInnes } 35952392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 360fb2e594dSBarry Smith EXTERN_C_BEGIN 3615615d1e5SSatish Balay #undef __FUNC__ 3625615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_TR" 363f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes ) 3644800dd8cSBarry Smith { 3656831982aSBarry Smith SNES_EQ_TR *neP; 3664800dd8cSBarry Smith 3673a40ed3dSBarry Smith PetscFunctionBegin; 3683a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 369a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 3703a40ed3dSBarry Smith } 371f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_TR; 372f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_TR; 373f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_TR; 374f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_TR; 375f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_TR; 376f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_TR; 377f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_TR; 3785baf8537SBarry Smith snes->nwork = 0; 379fbe28522SBarry Smith 3806831982aSBarry Smith neP = PetscNew(SNES_EQ_TR);CHKPTRQ(neP); 3816831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_TR)); 382fbe28522SBarry Smith snes->data = (void *) neP; 383fbe28522SBarry Smith neP->mu = 0.25; 384fbe28522SBarry Smith neP->eta = 0.75; 385fbe28522SBarry Smith neP->delta = 0.0; 386fbe28522SBarry Smith neP->delta0 = 0.2; 387fbe28522SBarry Smith neP->delta1 = 0.3; 388fbe28522SBarry Smith neP->delta2 = 0.75; 389fbe28522SBarry Smith neP->delta3 = 2.0; 390fbe28522SBarry Smith neP->sigma = 0.0001; 391fbe28522SBarry Smith neP->itflag = 0; 39252392280SLois Curfman McInnes neP->rnorm0 = 0; 39352392280SLois Curfman McInnes neP->ttol = 0; 3943a40ed3dSBarry Smith PetscFunctionReturn(0); 3954800dd8cSBarry Smith } 396fb2e594dSBarry Smith EXTERN_C_END 39782bf6240SBarry Smith 398