1*b2863d3aSBarry Smith /*$Id: ls.c,v 1.151 2000/02/02 20:10:06 bsmith Exp bsmith $*/ 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 604d965bbSLois Curfman McInnes 704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 804d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 1004d965bbSLois Curfman McInnes respectively. 1104d965bbSLois Curfman McInnes 12fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 1304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 1404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 15fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 1604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 1704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 1804d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 1904d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 2004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 2104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 2204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 2304d965bbSLois Curfman McInnes for all nonlinear solvers. 2404d965bbSLois Curfman McInnes 2504d965bbSLois Curfman McInnes Another key routine is: 2604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 2704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 2804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 2904d965bbSLois Curfman McInnes SNESSolve() if necessary. 3004d965bbSLois Curfman McInnes 3104d965bbSLois Curfman McInnes Additional basic routines are: 3204d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 3304d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 3404d965bbSLois Curfman McInnes have actually been used. 3504d965bbSLois Curfman McInnes These are called by application codes via the interface routines 3604d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 3704d965bbSLois Curfman McInnes 3804d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 3904d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 4004d965bbSLois Curfman McInnes above description applies to these categories also. 4104d965bbSLois Curfman McInnes 4204d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 435e42470aSBarry Smith /* 4404d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 4504d965bbSLois Curfman McInnes method with a line search. 465e76c431SBarry Smith 4704d965bbSLois Curfman McInnes Input Parameters: 4804d965bbSLois Curfman McInnes . snes - the SNES context 495e76c431SBarry Smith 5004d965bbSLois Curfman McInnes Output Parameter: 51c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 5204d965bbSLois Curfman McInnes 5304d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 545e76c431SBarry Smith 555e76c431SBarry Smith Notes: 565e76c431SBarry Smith This implements essentially a truncated Newton method with a 575e76c431SBarry Smith line search. By default a cubic backtracking line search 585e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 595e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 60393d2d9aSLois Curfman McInnes and Schnabel. 615e76c431SBarry Smith */ 625615d1e5SSatish Balay #undef __FUNC__ 63*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 64f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 655e76c431SBarry Smith { 666831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 6784c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 68112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 69329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 705e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 715e76c431SBarry Smith 723a40ed3dSBarry Smith PetscFunctionBegin; 73184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 74184914b5SBarry Smith 755e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 765e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 7739e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 785e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 795e42470aSBarry Smith G = snes->work[1]; 805e42470aSBarry Smith W = snes->work[2]; 815e76c431SBarry Smith 824c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 834c49b128SBarry Smith snes->iter = 0; 844c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 855334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 86cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 870f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 885e42470aSBarry Smith snes->norm = fnorm; 890f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 9084c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 9194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 925e76c431SBarry Smith 93184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 9494a9d846SBarry Smith 95d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 96d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 97d034289bSBarry Smith 985e76c431SBarry Smith for (i=0; i<maxits; i++) { 995e76c431SBarry Smith 100ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 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); 10378b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 1043505fcc1SBarry Smith snes->linear_its += lits; 105af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 106ea4d3ed3SLois Curfman McInnes 107ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 108ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 109ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 110ea4d3ed3SLois Curfman McInnes */ 11181b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 112af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 113af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 1143505fcc1SBarry Smith if (lsfail) { 1153505fcc1SBarry Smith snes->nfailures++; 1163505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 1173505fcc1SBarry Smith break; 1183505fcc1SBarry Smith } 1195e76c431SBarry Smith 12039e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 12139e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1225e76c431SBarry Smith fnorm = gnorm; 1235e76c431SBarry Smith 1240f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 12584c569c9SBarry Smith snes->iter = i+1; 1265e42470aSBarry Smith snes->norm = fnorm; 1270f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 12884c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 12994a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1305e76c431SBarry Smith 1315e76c431SBarry Smith /* Test for convergence */ 13229e0b56fSBarry Smith if (snes->converged) { 13329e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1343505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1353505fcc1SBarry Smith if (snes->reason) { 13616c95cacSBarry Smith break; 13716c95cacSBarry Smith } 13816c95cacSBarry Smith } 13929e0b56fSBarry Smith } 14039e2f89bSBarry Smith if (X != snes->vec_sol) { 141393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 14239e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 14339e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 14439e2f89bSBarry Smith } 14552392280SLois Curfman McInnes if (i == maxits) { 146981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 14752392280SLois Curfman McInnes i--; 1483505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 14952392280SLois Curfman McInnes } 1500f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1510f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1525e42470aSBarry Smith *outits = i+1; 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 1545e76c431SBarry Smith } 15504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 15604d965bbSLois Curfman McInnes /* 15704d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 1586831982aSBarry Smith of the SNESEQLS nonlinear solver. 15904d965bbSLois Curfman McInnes 16004d965bbSLois Curfman McInnes Input Parameter: 16104d965bbSLois Curfman McInnes . snes - the SNES context 16204d965bbSLois Curfman McInnes . x - the solution vector 16304d965bbSLois Curfman McInnes 16404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 16504d965bbSLois Curfman McInnes 16604d965bbSLois Curfman McInnes Notes: 16704d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 16804d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 16904d965bbSLois Curfman McInnes the call to SNESSolve(). 17004d965bbSLois Curfman McInnes */ 1715615d1e5SSatish Balay #undef __FUNC__ 172*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 173f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1745e76c431SBarry Smith { 1755e42470aSBarry Smith int ierr; 1763a40ed3dSBarry Smith 1773a40ed3dSBarry Smith PetscFunctionBegin; 17881b6cf68SLois Curfman McInnes snes->nwork = 4; 179d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1805e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 18181b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1823a40ed3dSBarry Smith PetscFunctionReturn(0); 1835e76c431SBarry Smith } 18404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 18504d965bbSLois Curfman McInnes /* 1866831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 18704d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 18804d965bbSLois Curfman McInnes 18904d965bbSLois Curfman McInnes Input Parameter: 19004d965bbSLois Curfman McInnes . snes - the SNES context 19104d965bbSLois Curfman McInnes 192de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 19304d965bbSLois Curfman McInnes */ 1945615d1e5SSatish Balay #undef __FUNC__ 195*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 196e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1975e76c431SBarry Smith { 198393d2d9aSLois Curfman McInnes int ierr; 1993a40ed3dSBarry Smith 2003a40ed3dSBarry Smith PetscFunctionBegin; 2015baf8537SBarry Smith if (snes->nwork) { 2024b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 20321c89e3eSBarry Smith } 2045c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2053a40ed3dSBarry Smith PetscFunctionReturn(0); 2065e76c431SBarry Smith } 20704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2085615d1e5SSatish Balay #undef __FUNC__ 209*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 21082bf6240SBarry Smith 2114b828684SBarry Smith /*@C 2125e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2135e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2145e76c431SBarry Smith to serve as a template and is not recommended for general use. 2155e76c431SBarry Smith 216c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 217c7afd0dbSLois Curfman McInnes 2185e76c431SBarry Smith Input Parameters: 219c7afd0dbSLois Curfman McInnes + snes - nonlinear context 220af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2215e76c431SBarry Smith . x - current iterate 2225e76c431SBarry Smith . f - residual evaluated at x 2235e76c431SBarry Smith . y - search direction (contains new iterate on output) 2245e76c431SBarry Smith . w - work vector 225c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2265e76c431SBarry Smith 227c4a48953SLois Curfman McInnes Output Parameters: 228c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2295e76c431SBarry Smith . y - new iterate (contains search direction on input) 2305e76c431SBarry Smith . gnorm - 2-norm of g 2315e76c431SBarry Smith . ynorm - 2-norm of search length 232c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 233fee21e36SBarry Smith 234c4a48953SLois Curfman McInnes Options Database Key: 235c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 236c4a48953SLois Curfman McInnes 23736851e7fSLois Curfman McInnes Level: advanced 23836851e7fSLois Curfman McInnes 23928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 24028ae5a21SLois Curfman McInnes 241f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 242af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2435e76c431SBarry Smith @*/ 244af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 245329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 2465e76c431SBarry Smith { 2475e42470aSBarry Smith int ierr; 2485334005bSBarry Smith Scalar mone = -1.0; 2496831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 2508f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2515334005bSBarry Smith 2523a40ed3dSBarry Smith PetscFunctionBegin; 253761aaf1bSLois Curfman McInnes *flag = 0; 2547857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 255a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 256ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 2578f99978dSLois Curfman McInnes if (neP->CheckStep) { 2588f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 2598f99978dSLois Curfman McInnes } 260ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 261a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 2627857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 2645e76c431SBarry Smith } 26504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2665615d1e5SSatish Balay #undef __FUNC__ 267*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 26882bf6240SBarry Smith 26929e0b56fSBarry Smith /*@C 27029e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 27129e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 27229e0b56fSBarry Smith even compute the norm of the function or search direction; this 27329e0b56fSBarry Smith is intended only when you know the full step is fine and are 27429e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 275c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 276c7afd0dbSLois Curfman McInnes 277c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 27829e0b56fSBarry Smith 27929e0b56fSBarry Smith Input Parameters: 280c7afd0dbSLois Curfman McInnes + snes - nonlinear context 281af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 28229e0b56fSBarry Smith . x - current iterate 28329e0b56fSBarry Smith . f - residual evaluated at x 28429e0b56fSBarry Smith . y - search direction (contains new iterate on output) 28529e0b56fSBarry Smith . w - work vector 286c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 28729e0b56fSBarry Smith 28829e0b56fSBarry Smith Output Parameters: 289c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 29029e0b56fSBarry Smith . gnorm - not changed 29129e0b56fSBarry Smith . ynorm - not changed 292c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 293fee21e36SBarry Smith 29429e0b56fSBarry Smith Options Database Key: 295c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 29629e0b56fSBarry Smith 2978cbba510SBarry Smith Notes: 298ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 299ea56f5baSLois Curfman McInnes either the options 300ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 301ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 302329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 303329f5518SBarry Smith otherwise, the SNES solver will generate an error. 304329f5518SBarry Smith 305329f5518SBarry Smith During the final iteration this will not evaluate the function at 306329f5518SBarry Smith the solution point. This is to save a function evaluation while 307329f5518SBarry Smith using pseudo-timestepping. 3088cbba510SBarry Smith 309ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 310ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 311ea56f5baSLois Curfman McInnes correct, since they are not computed. 312ea56f5baSLois Curfman McInnes 313ea56f5baSLois Curfman McInnes Level: advanced 3148cbba510SBarry Smith 31529e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 31629e0b56fSBarry Smith 31729e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 31829e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 31929e0b56fSBarry Smith @*/ 320af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 321329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 32229e0b56fSBarry Smith { 32329e0b56fSBarry Smith int ierr; 32429e0b56fSBarry Smith Scalar mone = -1.0; 3256831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3268f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 32729e0b56fSBarry Smith 3283a40ed3dSBarry Smith PetscFunctionBegin; 3298cbba510SBarry Smith *flag = 0; 33029e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 33129e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3328f99978dSLois Curfman McInnes if (neP->CheckStep) { 3338f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3348f99978dSLois Curfman McInnes } 335329f5518SBarry Smith 336329f5518SBarry Smith /* don't evaluate function the last time through */ 337329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 33829e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 339329f5518SBarry Smith } 34029e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3413a40ed3dSBarry Smith PetscFunctionReturn(0); 34229e0b56fSBarry Smith } 34304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34429e0b56fSBarry Smith #undef __FUNC__ 345*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 3464b828684SBarry Smith /*@C 347f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3485e76c431SBarry Smith 349c7afd0dbSLois Curfman McInnes Collective on SNES 350c7afd0dbSLois Curfman McInnes 3515e76c431SBarry Smith Input Parameters: 352c7afd0dbSLois Curfman McInnes + snes - nonlinear context 353af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3545e76c431SBarry Smith . x - current iterate 3555e76c431SBarry Smith . f - residual evaluated at x 3565e76c431SBarry Smith . y - search direction (contains new iterate on output) 3575e76c431SBarry Smith . w - work vector 358c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3595e76c431SBarry Smith 360393d2d9aSLois Curfman McInnes Output Parameters: 361c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3625e76c431SBarry Smith . y - new iterate (contains search direction on input) 3635e76c431SBarry Smith . gnorm - 2-norm of g 3645e76c431SBarry Smith . ynorm - 2-norm of search length 365c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 366fee21e36SBarry Smith 367c4a48953SLois Curfman McInnes Options Database Key: 368c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 369c4a48953SLois Curfman McInnes 3705e76c431SBarry Smith Notes: 3715e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3725e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3735e76c431SBarry Smith 37436851e7fSLois Curfman McInnes Level: advanced 37536851e7fSLois Curfman McInnes 37628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 37728ae5a21SLois Curfman McInnes 378af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3795e76c431SBarry Smith @*/ 380af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 381329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3825e76c431SBarry Smith { 383406556e6SLois Curfman McInnes /* 384406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 385406556e6SLois Curfman McInnes minimization problem: 386406556e6SLois Curfman McInnes min z(x): R^n -> R, 387406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 388406556e6SLois Curfman McInnes */ 389406556e6SLois Curfman McInnes 390329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 391329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 392aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 3935e42470aSBarry Smith Scalar cinitslope,clambda; 3946b5873e3SBarry Smith #endif 3955e42470aSBarry Smith int ierr,count; 3966831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3975334005bSBarry Smith Scalar mone = -1.0,scale; 3988f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3995e76c431SBarry Smith 4003a40ed3dSBarry Smith PetscFunctionBegin; 4017857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 402761aaf1bSLois Curfman McInnes *flag = 0; 4035e76c431SBarry Smith alpha = neP->alpha; 4045e76c431SBarry Smith maxstep = neP->maxstep; 4055e76c431SBarry Smith steptol = neP->steptol; 4065e76c431SBarry Smith 407cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 408a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 409a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 410a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 411a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 412a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 413ad922ac9SBarry Smith goto theend1; 41494a9d846SBarry Smith } 4155e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4165e42470aSBarry Smith scale = maxstep/(*ynorm); 417aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 418329f5518SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscRealPart(scale)); 4196b5873e3SBarry Smith #else 42094a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 4216b5873e3SBarry Smith #endif 422393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4235e76c431SBarry Smith *ynorm = maxstep; 4245e76c431SBarry Smith } 4255e76c431SBarry Smith minlambda = steptol/(*ynorm); 426a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 427aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 428a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 429329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 4305e42470aSBarry Smith #else 431a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4325e42470aSBarry Smith #endif 4335e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4345e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4355e76c431SBarry Smith 436393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4375334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 43878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 439313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 440406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 441393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 44294a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 443ad922ac9SBarry Smith goto theend1; 4445e76c431SBarry Smith } 4455e76c431SBarry Smith 4465e76c431SBarry Smith /* Fit points with quadratic */ 447313b5538SBarry Smith lambda = 1.0; 448a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4495e76c431SBarry Smith lambdaprev = lambda; 4505e76c431SBarry Smith gnormprev = *gnorm; 451329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 452ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 453ddd12767SBarry Smith else lambda = lambdatemp; 454393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 455ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 456aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 457ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4585e42470aSBarry Smith #else 459ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4605e42470aSBarry Smith #endif 46178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 462cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 463406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 464393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 465ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 466ad922ac9SBarry Smith goto theend1; 4675e76c431SBarry Smith } 4685e76c431SBarry Smith 4695e76c431SBarry Smith /* Fit points with cubic */ 4705e76c431SBarry Smith count = 1; 4715e76c431SBarry Smith while (1) { 4725e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4732b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4742b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 475ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 476393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 477761aaf1bSLois Curfman McInnes *flag = -1; break; 4785e76c431SBarry Smith } 479406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 480406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 481ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4822b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4835e76c431SBarry Smith d = b*b - 3*a*initslope; 4845e76c431SBarry Smith if (d < 0.0) d = 0.0; 4855e76c431SBarry Smith if (a == 0.0) { 4865e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4875e76c431SBarry Smith } else { 4885e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4895e76c431SBarry Smith } 4905e76c431SBarry Smith lambdaprev = lambda; 4915e76c431SBarry Smith gnormprev = *gnorm; 492329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 493329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 4945e76c431SBarry Smith else lambda = lambdatemp; 495393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 496ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 497aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 498ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 499393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5005e42470aSBarry Smith #else 501ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5025e42470aSBarry Smith #endif 50378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 504cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 505406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 506393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 507ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5082715565aSLois Curfman McInnes goto theend1; 5092b022350SLois Curfman McInnes } else { 5102b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5115e76c431SBarry Smith } 5125e76c431SBarry Smith count++; 5135e76c431SBarry Smith } 514ad922ac9SBarry Smith theend1: 5158f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5168f99978dSLois Curfman McInnes if (neP->CheckStep) { 5178f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5188f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5198f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5208f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5218f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5228f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5238f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5248f99978dSLois Curfman McInnes } 5258f99978dSLois Curfman McInnes } 5267857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5273a40ed3dSBarry Smith PetscFunctionReturn(0); 5285e76c431SBarry Smith } 52904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5305615d1e5SSatish Balay #undef __FUNC__ 531*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 5324b828684SBarry Smith /*@C 533f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5345e76c431SBarry Smith 535c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 536c7afd0dbSLois Curfman McInnes 5375e42470aSBarry Smith Input Parameters: 538c7afd0dbSLois Curfman McInnes + snes - the SNES context 539af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5405e42470aSBarry Smith . x - current iterate 5415e42470aSBarry Smith . f - residual evaluated at x 5425e42470aSBarry Smith . y - search direction (contains new iterate on output) 5435e42470aSBarry Smith . w - work vector 544c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5455e42470aSBarry Smith 546c4a48953SLois Curfman McInnes Output Parameters: 547c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5485e42470aSBarry Smith . y - new iterate (contains search direction on input) 5495e42470aSBarry Smith . gnorm - 2-norm of g 5505e42470aSBarry Smith . ynorm - 2-norm of search length 551c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 552fee21e36SBarry Smith 553c4a48953SLois Curfman McInnes Options Database Key: 554c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 555c4a48953SLois Curfman McInnes 5565e42470aSBarry Smith Notes: 5576831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 5585e42470aSBarry Smith 55936851e7fSLois Curfman McInnes Level: advanced 56036851e7fSLois Curfman McInnes 561f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 562f59ffdedSLois Curfman McInnes 563af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5645e42470aSBarry Smith @*/ 565af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 566329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 5675e76c431SBarry Smith { 568406556e6SLois Curfman McInnes /* 569406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 570406556e6SLois Curfman McInnes minimization problem: 571406556e6SLois Curfman McInnes min z(x): R^n -> R, 572406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 573406556e6SLois Curfman McInnes */ 574329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 575aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5765e42470aSBarry Smith Scalar cinitslope,clambda; 5776b5873e3SBarry Smith #endif 5785e42470aSBarry Smith int ierr,count; 5796831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 5805334005bSBarry Smith Scalar mone = -1.0,scale; 5818f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5825e76c431SBarry Smith 5833a40ed3dSBarry Smith PetscFunctionBegin; 5847857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 585761aaf1bSLois Curfman McInnes *flag = 0; 5865e42470aSBarry Smith alpha = neP->alpha; 5875e42470aSBarry Smith maxstep = neP->maxstep; 5885e42470aSBarry Smith steptol = neP->steptol; 5895e76c431SBarry Smith 5903505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 591b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 59294a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 593b37302c6SLois Curfman McInnes *gnorm = fnorm; 594b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 595b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 596ad922ac9SBarry Smith goto theend2; 59794a9d846SBarry Smith } 5985e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5995e42470aSBarry Smith scale = maxstep/(*ynorm); 600393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6015e42470aSBarry Smith *ynorm = maxstep; 6025e76c431SBarry Smith } 6035e42470aSBarry Smith minlambda = steptol/(*ynorm); 604a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 605aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 606a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 607329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6085e42470aSBarry Smith #else 609a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6105e42470aSBarry Smith #endif 6115e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6125e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6135e42470aSBarry Smith 614393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6155334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 61678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 617cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 618406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 619393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 62094a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 621ad922ac9SBarry Smith goto theend2; 6225e42470aSBarry Smith } 6235e42470aSBarry Smith 6245e42470aSBarry Smith /* Fit points with quadratic */ 625313b5538SBarry Smith lambda = 1.0; 6265e42470aSBarry Smith count = 1; 6275e42470aSBarry Smith while (1) { 6285e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 629981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 630981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 631ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 632393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 633761aaf1bSLois Curfman McInnes *flag = -1; break; 6345e42470aSBarry Smith } 635a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 636329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 637329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 638329f5518SBarry Smith else lambda = lambdatemp; 639393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6403505fcc1SBarry Smith lambdaneg = -lambda; 641aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6423505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6435e42470aSBarry Smith #else 6443505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6455e42470aSBarry Smith #endif 64678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 647cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 648406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 649393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 650981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 65106259719SBarry Smith break; 6525e42470aSBarry Smith } 6535e42470aSBarry Smith count++; 6545e42470aSBarry Smith } 655ad922ac9SBarry Smith theend2: 6568f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6578f99978dSLois Curfman McInnes if (neP->CheckStep) { 6588f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6598f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6608f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6618f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6628f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6638f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6648f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6658f99978dSLois Curfman McInnes } 6668f99978dSLois Curfman McInnes } 6677857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6683a40ed3dSBarry Smith PetscFunctionReturn(0); 6695e76c431SBarry Smith } 67004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6715615d1e5SSatish Balay #undef __FUNC__ 672*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 673c9e6a524SLois Curfman McInnes /*@C 67477c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 6756831982aSBarry Smith by the method SNESEQLS. 6765e76c431SBarry Smith 6775e76c431SBarry Smith Input Parameters: 678c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 679af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 680c7afd0dbSLois Curfman McInnes - func - pointer to int function 6815e76c431SBarry Smith 682fee21e36SBarry Smith Collective on SNES 683fee21e36SBarry Smith 684c4a48953SLois Curfman McInnes Available Routines: 685c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 686f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 687f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 688af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6895e76c431SBarry Smith 690c4a48953SLois Curfman McInnes Options Database Keys: 691af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 692c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 693c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 694c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 695c4a48953SLois Curfman McInnes 6965e76c431SBarry Smith Calling sequence of func: 697bd208895SLois Curfman McInnes .vb 698af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 699329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 700bd208895SLois Curfman McInnes .ve 7015e76c431SBarry Smith 7025e76c431SBarry Smith Input parameters for func: 703c7afd0dbSLois Curfman McInnes + snes - nonlinear context 704af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7055e76c431SBarry Smith . x - current iterate 7065e76c431SBarry Smith . f - residual evaluated at x 7075e76c431SBarry Smith . y - search direction (contains new iterate on output) 7085e76c431SBarry Smith . w - work vector 709c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7105e76c431SBarry Smith 7115e76c431SBarry Smith Output parameters for func: 712c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7135e76c431SBarry Smith . y - new iterate (contains search direction on input) 7145e76c431SBarry Smith . gnorm - 2-norm of g 7155e76c431SBarry Smith . ynorm - 2-norm of search length 716c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 717761aaf1bSLois Curfman McInnes on failure. 718f59ffdedSLois Curfman McInnes 71936851e7fSLois Curfman McInnes Level: advanced 72036851e7fSLois Curfman McInnes 721f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 722f59ffdedSLois Curfman McInnes 7234ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7245e76c431SBarry Smith @*/ 725329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7265e76c431SBarry Smith { 727329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 72882bf6240SBarry Smith 7293a40ed3dSBarry Smith PetscFunctionBegin; 73094d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 73182bf6240SBarry Smith if (f) { 732af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 73382bf6240SBarry Smith } 7343a40ed3dSBarry Smith PetscFunctionReturn(0); 7355e76c431SBarry Smith } 73604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 737fb2e594dSBarry Smith EXTERN_C_BEGIN 73882bf6240SBarry Smith #undef __FUNC__ 739*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 740af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 741af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 74282bf6240SBarry Smith { 74382bf6240SBarry Smith PetscFunctionBegin; 7446831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7456831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 74682bf6240SBarry Smith PetscFunctionReturn(0); 74782bf6240SBarry Smith } 748fb2e594dSBarry Smith EXTERN_C_END 74904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 750c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 751*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 752c8dd0c0dSLois Curfman McInnes /*@C 753530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7546831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 755c8dd0c0dSLois Curfman McInnes 756c8dd0c0dSLois Curfman McInnes Input Parameters: 757c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 758c8dd0c0dSLois Curfman McInnes . func - pointer to int function 759c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 760c8dd0c0dSLois Curfman McInnes 761c8dd0c0dSLois Curfman McInnes Collective on SNES 762c8dd0c0dSLois Curfman McInnes 763c8dd0c0dSLois Curfman McInnes Calling sequence of func: 764c8dd0c0dSLois Curfman McInnes .vb 765b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 766c8dd0c0dSLois Curfman McInnes .ve 767b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 768b1ae27eaSLois Curfman McInnes on failure. 769c8dd0c0dSLois Curfman McInnes 770c8dd0c0dSLois Curfman McInnes Input parameters for func: 771c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 772c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7738f99978dSLois Curfman McInnes - x - current candidate iterate 774c8dd0c0dSLois Curfman McInnes 775c8dd0c0dSLois Curfman McInnes Output parameters for func: 776c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 777c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 778c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 779c8dd0c0dSLois Curfman McInnes 780c8dd0c0dSLois Curfman McInnes Level: advanced 781c8dd0c0dSLois Curfman McInnes 7828f99978dSLois Curfman McInnes Notes: 783b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 784b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 785b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 786b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 787ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 7888f99978dSLois Curfman McInnes 789b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 790b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 791b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 792ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 793ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 794ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 795ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 796b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 7978f99978dSLois Curfman McInnes 798c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 799c8dd0c0dSLois Curfman McInnes 800c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 801c8dd0c0dSLois Curfman McInnes @*/ 8028f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 803c8dd0c0dSLois Curfman McInnes { 8048f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 805c8dd0c0dSLois Curfman McInnes 806c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 807c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 808c8dd0c0dSLois Curfman McInnes if (f) { 809c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 810c8dd0c0dSLois Curfman McInnes } 811c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 812c8dd0c0dSLois Curfman McInnes } 813c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 814c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 815c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 816*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8178f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 818c8dd0c0dSLois Curfman McInnes { 819c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8206831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8216831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 822c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 823c8dd0c0dSLois Curfman McInnes } 824c8dd0c0dSLois Curfman McInnes EXTERN_C_END 825c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 82604d965bbSLois Curfman McInnes /* 8276831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 82882bf6240SBarry Smith 82904d965bbSLois Curfman McInnes Input Parameter: 83004d965bbSLois Curfman McInnes . snes - the SNES context 83104d965bbSLois Curfman McInnes 83204d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 83304d965bbSLois Curfman McInnes */ 8345615d1e5SSatish Balay #undef __FUNC__ 835*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 836f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8375e42470aSBarry Smith { 8386831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 839d132466eSBarry Smith int ierr; 8406daaf66cSBarry Smith 8413a40ed3dSBarry Smith PetscFunctionBegin; 8426831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 843d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 844d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 845d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 846d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8473a40ed3dSBarry Smith PetscFunctionReturn(0); 8485e42470aSBarry Smith } 84904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 85004d965bbSLois Curfman McInnes /* 8516831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 85204d965bbSLois Curfman McInnes 85304d965bbSLois Curfman McInnes Input Parameters: 85404d965bbSLois Curfman McInnes . SNES - the SNES context 85504d965bbSLois Curfman McInnes . viewer - visualization context 85604d965bbSLois Curfman McInnes 85704d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 85804d965bbSLois Curfman McInnes */ 8595615d1e5SSatish Balay #undef __FUNC__ 860*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 861e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 862a935fc98SLois Curfman McInnes { 8636831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 86419bcc07fSBarry Smith char *cstr; 86551695f54SLois Curfman McInnes int ierr; 8666831982aSBarry Smith PetscTruth isascii; 867a935fc98SLois Curfman McInnes 8683a40ed3dSBarry Smith PetscFunctionBegin; 8696831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8700f5bd95cSBarry Smith if (isascii) { 87119bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 87219bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 87319bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 87419bcc07fSBarry Smith else cstr = "unknown"; 875d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 876d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 8775cd90555SBarry Smith } else { 8780f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 87919bcc07fSBarry Smith } 8803a40ed3dSBarry Smith PetscFunctionReturn(0); 881a935fc98SLois Curfman McInnes } 88204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 88304d965bbSLois Curfman McInnes /* 8846831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 88504d965bbSLois Curfman McInnes 88604d965bbSLois Curfman McInnes Input Parameter: 88704d965bbSLois Curfman McInnes . snes - the SNES context 88804d965bbSLois Curfman McInnes 88904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 89004d965bbSLois Curfman McInnes */ 8915615d1e5SSatish Balay #undef __FUNC__ 892*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 893f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8945e42470aSBarry Smith { 8956831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 8965e42470aSBarry Smith char ver[16]; 8975e42470aSBarry Smith double tmp; 898f1af5d2fSBarry Smith int ierr; 899f1af5d2fSBarry Smith PetscTruth flg; 9005e42470aSBarry Smith 9013a40ed3dSBarry Smith PetscFunctionBegin; 90209d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 90325cdf11fSBarry Smith if (flg) { 9045e42470aSBarry Smith ls->alpha = tmp; 9055e42470aSBarry Smith } 90609d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 90725cdf11fSBarry Smith if (flg) { 9085e42470aSBarry Smith ls->maxstep = tmp; 9095e42470aSBarry Smith } 91009d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 91125cdf11fSBarry Smith if (flg) { 9125e42470aSBarry Smith ls->steptol = tmp; 9135e42470aSBarry Smith } 91409d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 91525cdf11fSBarry Smith if (flg) { 916c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9170f5bd95cSBarry Smith 918c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 919c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 920c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 921c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9220f5bd95cSBarry Smith 9230f5bd95cSBarry Smith if (isbasic) { 924af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9250f5bd95cSBarry Smith } else if (isnonorms) { 926af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9270f5bd95cSBarry Smith } else if (isquad) { 928af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9290f5bd95cSBarry Smith } else if (iscubic) { 930af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9315e42470aSBarry Smith } 932a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9335e42470aSBarry Smith } 9343a40ed3dSBarry Smith PetscFunctionReturn(0); 9355e42470aSBarry Smith } 93604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 93704d965bbSLois Curfman McInnes /* 9386831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9396831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 94004d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 94104d965bbSLois Curfman McInnes 94204d965bbSLois Curfman McInnes 94304d965bbSLois Curfman McInnes Input Parameter: 94404d965bbSLois Curfman McInnes . snes - the SNES context 94504d965bbSLois Curfman McInnes 94604d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 94704d965bbSLois Curfman McInnes */ 948fb2e594dSBarry Smith EXTERN_C_BEGIN 9495615d1e5SSatish Balay #undef __FUNC__ 950*b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 951f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9525e42470aSBarry Smith { 95382bf6240SBarry Smith int ierr; 9546831982aSBarry Smith SNES_EQ_LS *neP; 9555e42470aSBarry Smith 9563a40ed3dSBarry Smith PetscFunctionBegin; 957a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 958a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 959a8c6a408SBarry Smith } 96082bf6240SBarry Smith 961f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 962f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 963f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 964f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 965f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 966f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 967f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9685baf8537SBarry Smith snes->nwork = 0; 9695e42470aSBarry Smith 9706831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 9716831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 9725e42470aSBarry Smith snes->data = (void*)neP; 9735e42470aSBarry Smith neP->alpha = 1.e-4; 9745e42470aSBarry Smith neP->maxstep = 1.e8; 9755e42470aSBarry Smith neP->steptol = 1.e-12; 9765e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 977c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 978c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 979c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 98082bf6240SBarry Smith 981f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 982e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 983f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 984c8dd0c0dSLois Curfman McInnes (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 98582bf6240SBarry Smith 9863a40ed3dSBarry Smith PetscFunctionReturn(0); 9875e42470aSBarry Smith } 988fb2e594dSBarry Smith EXTERN_C_END 98982bf6240SBarry Smith 99082bf6240SBarry Smith 99182bf6240SBarry Smith 992