1*0c97f09dSSatish Balay /*$Id: ls.c,v 1.154 2000/04/18 23:03:49 bsmith Exp balay $*/ 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__ 63b2863d3aSBarry 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); 142e15f7bb5SBarry Smith } 143e15f7bb5SBarry Smith if (F != snes->vec_func) { 144e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 145e15f7bb5SBarry Smith } 14639e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 14739e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 14852392280SLois Curfman McInnes if (i == maxits) { 149981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 15052392280SLois Curfman McInnes i--; 1513505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 15252392280SLois Curfman McInnes } 1530f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1540f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1555e42470aSBarry Smith *outits = i+1; 1563a40ed3dSBarry Smith PetscFunctionReturn(0); 1575e76c431SBarry Smith } 15804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 15904d965bbSLois Curfman McInnes /* 16004d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 1616831982aSBarry Smith of the SNESEQLS nonlinear solver. 16204d965bbSLois Curfman McInnes 16304d965bbSLois Curfman McInnes Input Parameter: 16404d965bbSLois Curfman McInnes . snes - the SNES context 16504d965bbSLois Curfman McInnes . x - the solution vector 16604d965bbSLois Curfman McInnes 16704d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 16804d965bbSLois Curfman McInnes 16904d965bbSLois Curfman McInnes Notes: 17004d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 17104d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 17204d965bbSLois Curfman McInnes the call to SNESSolve(). 17304d965bbSLois Curfman McInnes */ 1745615d1e5SSatish Balay #undef __FUNC__ 175b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 176f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1775e76c431SBarry Smith { 1785e42470aSBarry Smith int ierr; 1793a40ed3dSBarry Smith 1803a40ed3dSBarry Smith PetscFunctionBegin; 18181b6cf68SLois Curfman McInnes snes->nwork = 4; 182d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1835e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 18481b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1853a40ed3dSBarry Smith PetscFunctionReturn(0); 1865e76c431SBarry Smith } 18704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 18804d965bbSLois Curfman McInnes /* 1896831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 19004d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 19104d965bbSLois Curfman McInnes 19204d965bbSLois Curfman McInnes Input Parameter: 19304d965bbSLois Curfman McInnes . snes - the SNES context 19404d965bbSLois Curfman McInnes 195de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 19604d965bbSLois Curfman McInnes */ 1975615d1e5SSatish Balay #undef __FUNC__ 198b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 199e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2005e76c431SBarry Smith { 201393d2d9aSLois Curfman McInnes int ierr; 2023a40ed3dSBarry Smith 2033a40ed3dSBarry Smith PetscFunctionBegin; 2045baf8537SBarry Smith if (snes->nwork) { 2054b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 20621c89e3eSBarry Smith } 2075c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 2095e76c431SBarry Smith } 21004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2115615d1e5SSatish Balay #undef __FUNC__ 212b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 21382bf6240SBarry Smith 2144b828684SBarry Smith /*@C 2155e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2165e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2175e76c431SBarry Smith to serve as a template and is not recommended for general use. 2185e76c431SBarry Smith 219c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 220c7afd0dbSLois Curfman McInnes 2215e76c431SBarry Smith Input Parameters: 222c7afd0dbSLois Curfman McInnes + snes - nonlinear context 223af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2245e76c431SBarry Smith . x - current iterate 2255e76c431SBarry Smith . f - residual evaluated at x 2265e76c431SBarry Smith . y - search direction (contains new iterate on output) 2275e76c431SBarry Smith . w - work vector 228c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2295e76c431SBarry Smith 230c4a48953SLois Curfman McInnes Output Parameters: 231c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2325e76c431SBarry Smith . y - new iterate (contains search direction on input) 2335e76c431SBarry Smith . gnorm - 2-norm of g 2345e76c431SBarry Smith . ynorm - 2-norm of search length 235c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 236fee21e36SBarry Smith 237c4a48953SLois Curfman McInnes Options Database Key: 238c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 239c4a48953SLois Curfman McInnes 24036851e7fSLois Curfman McInnes Level: advanced 24136851e7fSLois Curfman McInnes 24228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 24328ae5a21SLois Curfman McInnes 244f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 245af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2465e76c431SBarry Smith @*/ 247af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 248329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 2495e76c431SBarry Smith { 2505e42470aSBarry Smith int ierr; 2515334005bSBarry Smith Scalar mone = -1.0; 2526831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 2538f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2545334005bSBarry Smith 2553a40ed3dSBarry Smith PetscFunctionBegin; 256761aaf1bSLois Curfman McInnes *flag = 0; 2577857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 258a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 259ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 2608f99978dSLois Curfman McInnes if (neP->CheckStep) { 2618f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 2628f99978dSLois Curfman McInnes } 263ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 264a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 2657857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2663a40ed3dSBarry Smith PetscFunctionReturn(0); 2675e76c431SBarry Smith } 26804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2695615d1e5SSatish Balay #undef __FUNC__ 270b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 27182bf6240SBarry Smith 27229e0b56fSBarry Smith /*@C 27329e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 27429e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 27529e0b56fSBarry Smith even compute the norm of the function or search direction; this 27629e0b56fSBarry Smith is intended only when you know the full step is fine and are 27729e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 278c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 279c7afd0dbSLois Curfman McInnes 280c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 28129e0b56fSBarry Smith 28229e0b56fSBarry Smith Input Parameters: 283c7afd0dbSLois Curfman McInnes + snes - nonlinear context 284af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 28529e0b56fSBarry Smith . x - current iterate 28629e0b56fSBarry Smith . f - residual evaluated at x 28729e0b56fSBarry Smith . y - search direction (contains new iterate on output) 28829e0b56fSBarry Smith . w - work vector 289c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 29029e0b56fSBarry Smith 29129e0b56fSBarry Smith Output Parameters: 292c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 29329e0b56fSBarry Smith . gnorm - not changed 29429e0b56fSBarry Smith . ynorm - not changed 295c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 296fee21e36SBarry Smith 29729e0b56fSBarry Smith Options Database Key: 298c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 29929e0b56fSBarry Smith 3008cbba510SBarry Smith Notes: 301ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 302ea56f5baSLois Curfman McInnes either the options 303ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 304ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 305329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 306329f5518SBarry Smith otherwise, the SNES solver will generate an error. 307329f5518SBarry Smith 308329f5518SBarry Smith During the final iteration this will not evaluate the function at 309329f5518SBarry Smith the solution point. This is to save a function evaluation while 310329f5518SBarry Smith using pseudo-timestepping. 3118cbba510SBarry Smith 312ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 313ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 314ea56f5baSLois Curfman McInnes correct, since they are not computed. 315ea56f5baSLois Curfman McInnes 316ea56f5baSLois Curfman McInnes Level: advanced 3178cbba510SBarry Smith 31829e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 31929e0b56fSBarry Smith 32029e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 32129e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 32229e0b56fSBarry Smith @*/ 323af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 324329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 32529e0b56fSBarry Smith { 32629e0b56fSBarry Smith int ierr; 32729e0b56fSBarry Smith Scalar mone = -1.0; 3286831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3298f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 33029e0b56fSBarry Smith 3313a40ed3dSBarry Smith PetscFunctionBegin; 3328cbba510SBarry Smith *flag = 0; 33329e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 33429e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3358f99978dSLois Curfman McInnes if (neP->CheckStep) { 3368f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3378f99978dSLois Curfman McInnes } 338329f5518SBarry Smith 339329f5518SBarry Smith /* don't evaluate function the last time through */ 340329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 34129e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 342329f5518SBarry Smith } 34329e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3443a40ed3dSBarry Smith PetscFunctionReturn(0); 34529e0b56fSBarry Smith } 34604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34729e0b56fSBarry Smith #undef __FUNC__ 348b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 3494b828684SBarry Smith /*@C 350f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3515e76c431SBarry Smith 352c7afd0dbSLois Curfman McInnes Collective on SNES 353c7afd0dbSLois Curfman McInnes 3545e76c431SBarry Smith Input Parameters: 355c7afd0dbSLois Curfman McInnes + snes - nonlinear context 356af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3575e76c431SBarry Smith . x - current iterate 3585e76c431SBarry Smith . f - residual evaluated at x 3595e76c431SBarry Smith . y - search direction (contains new iterate on output) 3605e76c431SBarry Smith . w - work vector 361c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3625e76c431SBarry Smith 363393d2d9aSLois Curfman McInnes Output Parameters: 364c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3655e76c431SBarry Smith . y - new iterate (contains search direction on input) 3665e76c431SBarry Smith . gnorm - 2-norm of g 3675e76c431SBarry Smith . ynorm - 2-norm of search length 368c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 369fee21e36SBarry Smith 370c4a48953SLois Curfman McInnes Options Database Key: 371c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 372c4a48953SLois Curfman McInnes 3735e76c431SBarry Smith Notes: 3745e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3755e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3765e76c431SBarry Smith 37736851e7fSLois Curfman McInnes Level: advanced 37836851e7fSLois Curfman McInnes 37928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 38028ae5a21SLois Curfman McInnes 381af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3825e76c431SBarry Smith @*/ 383af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 384329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3855e76c431SBarry Smith { 386406556e6SLois Curfman McInnes /* 387406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 388406556e6SLois Curfman McInnes minimization problem: 389406556e6SLois Curfman McInnes min z(x): R^n -> R, 390406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 391406556e6SLois Curfman McInnes */ 392406556e6SLois Curfman McInnes 393329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 394329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 395aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 3965e42470aSBarry Smith Scalar cinitslope,clambda; 3976b5873e3SBarry Smith #endif 3985e42470aSBarry Smith int ierr,count; 3996831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4005334005bSBarry Smith Scalar mone = -1.0,scale; 4018f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4025e76c431SBarry Smith 4033a40ed3dSBarry Smith PetscFunctionBegin; 4047857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 405761aaf1bSLois Curfman McInnes *flag = 0; 4065e76c431SBarry Smith alpha = neP->alpha; 4075e76c431SBarry Smith maxstep = neP->maxstep; 4085e76c431SBarry Smith steptol = neP->steptol; 4095e76c431SBarry Smith 410cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 411a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 412a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 413a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 414a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 415a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 416ad922ac9SBarry Smith goto theend1; 41794a9d846SBarry Smith } 4185e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4195e42470aSBarry Smith scale = maxstep/(*ynorm); 420aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 421329f5518SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscRealPart(scale)); 4226b5873e3SBarry Smith #else 42394a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 4246b5873e3SBarry Smith #endif 425393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4265e76c431SBarry Smith *ynorm = maxstep; 4275e76c431SBarry Smith } 4285e76c431SBarry Smith minlambda = steptol/(*ynorm); 429a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 430aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 431a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 432329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 4335e42470aSBarry Smith #else 434a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4355e42470aSBarry Smith #endif 4365e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4375e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4385e76c431SBarry Smith 439393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4405334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 44178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 442313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 443406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 444393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 44594a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 446ad922ac9SBarry Smith goto theend1; 4475e76c431SBarry Smith } 4485e76c431SBarry Smith 4495e76c431SBarry Smith /* Fit points with quadratic */ 450313b5538SBarry Smith lambda = 1.0; 451a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4525e76c431SBarry Smith lambdaprev = lambda; 4535e76c431SBarry Smith gnormprev = *gnorm; 454329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 455ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 456ddd12767SBarry Smith else lambda = lambdatemp; 457393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 458ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 459aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 460ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4615e42470aSBarry Smith #else 462ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4635e42470aSBarry Smith #endif 46478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 465cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 466406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 467393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 468ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 469ad922ac9SBarry Smith goto theend1; 4705e76c431SBarry Smith } 4715e76c431SBarry Smith 4725e76c431SBarry Smith /* Fit points with cubic */ 4735e76c431SBarry Smith count = 1; 4745e76c431SBarry Smith while (1) { 4755e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4762b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4772b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 478ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 479393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 480761aaf1bSLois Curfman McInnes *flag = -1; break; 4815e76c431SBarry Smith } 482406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 483406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 484ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4852b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4865e76c431SBarry Smith d = b*b - 3*a*initslope; 4875e76c431SBarry Smith if (d < 0.0) d = 0.0; 4885e76c431SBarry Smith if (a == 0.0) { 4895e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4905e76c431SBarry Smith } else { 4915e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4925e76c431SBarry Smith } 4935e76c431SBarry Smith lambdaprev = lambda; 4945e76c431SBarry Smith gnormprev = *gnorm; 495329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 496329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 4975e76c431SBarry Smith else lambda = lambdatemp; 498393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 499ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 500aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 501ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 502393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5035e42470aSBarry Smith #else 504ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5055e42470aSBarry Smith #endif 50678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 507cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 508406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 509393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 510ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5112715565aSLois Curfman McInnes goto theend1; 5122b022350SLois Curfman McInnes } else { 5132b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5145e76c431SBarry Smith } 5155e76c431SBarry Smith count++; 5165e76c431SBarry Smith } 517ad922ac9SBarry Smith theend1: 5188f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5198f99978dSLois Curfman McInnes if (neP->CheckStep) { 5208f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5218f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5228f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5238f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5248f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5258f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5268f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5278f99978dSLois Curfman McInnes } 5288f99978dSLois Curfman McInnes } 5297857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5303a40ed3dSBarry Smith PetscFunctionReturn(0); 5315e76c431SBarry Smith } 53204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5335615d1e5SSatish Balay #undef __FUNC__ 534b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 5354b828684SBarry Smith /*@C 536f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5375e76c431SBarry Smith 538c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 539c7afd0dbSLois Curfman McInnes 5405e42470aSBarry Smith Input Parameters: 541c7afd0dbSLois Curfman McInnes + snes - the SNES context 542af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5435e42470aSBarry Smith . x - current iterate 5445e42470aSBarry Smith . f - residual evaluated at x 5455e42470aSBarry Smith . y - search direction (contains new iterate on output) 5465e42470aSBarry Smith . w - work vector 547c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5485e42470aSBarry Smith 549c4a48953SLois Curfman McInnes Output Parameters: 550c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5515e42470aSBarry Smith . y - new iterate (contains search direction on input) 5525e42470aSBarry Smith . gnorm - 2-norm of g 5535e42470aSBarry Smith . ynorm - 2-norm of search length 554c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 555fee21e36SBarry Smith 556c4a48953SLois Curfman McInnes Options Database Key: 557c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 558c4a48953SLois Curfman McInnes 5595e42470aSBarry Smith Notes: 5606831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 5615e42470aSBarry Smith 56236851e7fSLois Curfman McInnes Level: advanced 56336851e7fSLois Curfman McInnes 564f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 565f59ffdedSLois Curfman McInnes 566af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5675e42470aSBarry Smith @*/ 568af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 569329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 5705e76c431SBarry Smith { 571406556e6SLois Curfman McInnes /* 572406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 573406556e6SLois Curfman McInnes minimization problem: 574406556e6SLois Curfman McInnes min z(x): R^n -> R, 575406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 576406556e6SLois Curfman McInnes */ 577329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 578aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5795e42470aSBarry Smith Scalar cinitslope,clambda; 5806b5873e3SBarry Smith #endif 5815e42470aSBarry Smith int ierr,count; 5826831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 5835334005bSBarry Smith Scalar mone = -1.0,scale; 5848f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5855e76c431SBarry Smith 5863a40ed3dSBarry Smith PetscFunctionBegin; 5877857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 588761aaf1bSLois Curfman McInnes *flag = 0; 5895e42470aSBarry Smith alpha = neP->alpha; 5905e42470aSBarry Smith maxstep = neP->maxstep; 5915e42470aSBarry Smith steptol = neP->steptol; 5925e76c431SBarry Smith 5933505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 594b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 59594a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 596b37302c6SLois Curfman McInnes *gnorm = fnorm; 597b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 598b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 599ad922ac9SBarry Smith goto theend2; 60094a9d846SBarry Smith } 6015e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6025e42470aSBarry Smith scale = maxstep/(*ynorm); 603393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6045e42470aSBarry Smith *ynorm = maxstep; 6055e76c431SBarry Smith } 6065e42470aSBarry Smith minlambda = steptol/(*ynorm); 607a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 608aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 609a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 610329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6115e42470aSBarry Smith #else 612a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6135e42470aSBarry Smith #endif 6145e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6155e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6165e42470aSBarry Smith 617393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6185334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 61978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 620cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 621406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 622393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 62394a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 624ad922ac9SBarry Smith goto theend2; 6255e42470aSBarry Smith } 6265e42470aSBarry Smith 6275e42470aSBarry Smith /* Fit points with quadratic */ 628313b5538SBarry Smith lambda = 1.0; 6295e42470aSBarry Smith count = 1; 6305e42470aSBarry Smith while (1) { 6315e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 632981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 633981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 634ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 635393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 636761aaf1bSLois Curfman McInnes *flag = -1; break; 6375e42470aSBarry Smith } 638a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 639329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 640329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 641329f5518SBarry Smith else lambda = lambdatemp; 642393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6433505fcc1SBarry Smith lambdaneg = -lambda; 644aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6453505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6465e42470aSBarry Smith #else 6473505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6485e42470aSBarry Smith #endif 64978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 650cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 651406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 652393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 653981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 65406259719SBarry Smith break; 6555e42470aSBarry Smith } 6565e42470aSBarry Smith count++; 6575e42470aSBarry Smith } 658ad922ac9SBarry Smith theend2: 6598f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6608f99978dSLois Curfman McInnes if (neP->CheckStep) { 6618f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6628f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6638f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6648f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6658f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6668f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6678f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6688f99978dSLois Curfman McInnes } 6698f99978dSLois Curfman McInnes } 6707857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6713a40ed3dSBarry Smith PetscFunctionReturn(0); 6725e76c431SBarry Smith } 67304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6745615d1e5SSatish Balay #undef __FUNC__ 675b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 676c9e6a524SLois Curfman McInnes /*@C 67777c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 6786831982aSBarry Smith by the method SNESEQLS. 6795e76c431SBarry Smith 6805e76c431SBarry Smith Input Parameters: 681c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 682af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 683c7afd0dbSLois Curfman McInnes - func - pointer to int function 6845e76c431SBarry Smith 685fee21e36SBarry Smith Collective on SNES 686fee21e36SBarry Smith 687c4a48953SLois Curfman McInnes Available Routines: 688c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 689f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 690f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 691af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6925e76c431SBarry Smith 693c4a48953SLois Curfman McInnes Options Database Keys: 694af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 695c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 696c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 697c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 698c4a48953SLois Curfman McInnes 6995e76c431SBarry Smith Calling sequence of func: 700bd208895SLois Curfman McInnes .vb 701af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 702329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 703bd208895SLois Curfman McInnes .ve 7045e76c431SBarry Smith 7055e76c431SBarry Smith Input parameters for func: 706c7afd0dbSLois Curfman McInnes + snes - nonlinear context 707af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7085e76c431SBarry Smith . x - current iterate 7095e76c431SBarry Smith . f - residual evaluated at x 7105e76c431SBarry Smith . y - search direction (contains new iterate on output) 7115e76c431SBarry Smith . w - work vector 712c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7135e76c431SBarry Smith 7145e76c431SBarry Smith Output parameters for func: 715c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7165e76c431SBarry Smith . y - new iterate (contains search direction on input) 7175e76c431SBarry Smith . gnorm - 2-norm of g 7185e76c431SBarry Smith . ynorm - 2-norm of search length 719c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 720761aaf1bSLois Curfman McInnes on failure. 721f59ffdedSLois Curfman McInnes 72236851e7fSLois Curfman McInnes Level: advanced 72336851e7fSLois Curfman McInnes 724f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 725f59ffdedSLois Curfman McInnes 7264ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7275e76c431SBarry Smith @*/ 728329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7295e76c431SBarry Smith { 730329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 73182bf6240SBarry Smith 7323a40ed3dSBarry Smith PetscFunctionBegin; 73394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 73482bf6240SBarry Smith if (f) { 735af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 73682bf6240SBarry Smith } 7373a40ed3dSBarry Smith PetscFunctionReturn(0); 7385e76c431SBarry Smith } 73904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 740fb2e594dSBarry Smith EXTERN_C_BEGIN 74182bf6240SBarry Smith #undef __FUNC__ 742b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 743af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 744af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 74582bf6240SBarry Smith { 74682bf6240SBarry Smith PetscFunctionBegin; 7476831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7486831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 74982bf6240SBarry Smith PetscFunctionReturn(0); 75082bf6240SBarry Smith } 751fb2e594dSBarry Smith EXTERN_C_END 75204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 753c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 754b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 755c8dd0c0dSLois Curfman McInnes /*@C 756530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7576831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 758c8dd0c0dSLois Curfman McInnes 759c8dd0c0dSLois Curfman McInnes Input Parameters: 760c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 761c8dd0c0dSLois Curfman McInnes . func - pointer to int function 762c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 763c8dd0c0dSLois Curfman McInnes 764c8dd0c0dSLois Curfman McInnes Collective on SNES 765c8dd0c0dSLois Curfman McInnes 766c8dd0c0dSLois Curfman McInnes Calling sequence of func: 767c8dd0c0dSLois Curfman McInnes .vb 768b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 769c8dd0c0dSLois Curfman McInnes .ve 770b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 771b1ae27eaSLois Curfman McInnes on failure. 772c8dd0c0dSLois Curfman McInnes 773c8dd0c0dSLois Curfman McInnes Input parameters for func: 774c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 775c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7768f99978dSLois Curfman McInnes - x - current candidate iterate 777c8dd0c0dSLois Curfman McInnes 778c8dd0c0dSLois Curfman McInnes Output parameters for func: 779c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 780c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 781c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 782c8dd0c0dSLois Curfman McInnes 783c8dd0c0dSLois Curfman McInnes Level: advanced 784c8dd0c0dSLois Curfman McInnes 7858f99978dSLois Curfman McInnes Notes: 786b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 787b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 788b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 789b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 790ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 7918f99978dSLois Curfman McInnes 792b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 793b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 794b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 795ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 796ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 797ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 798ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 799b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8008f99978dSLois Curfman McInnes 801c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 802c8dd0c0dSLois Curfman McInnes 803c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 804c8dd0c0dSLois Curfman McInnes @*/ 8058f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 806c8dd0c0dSLois Curfman McInnes { 8078f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 808c8dd0c0dSLois Curfman McInnes 809c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 810c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 811c8dd0c0dSLois Curfman McInnes if (f) { 812c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 813c8dd0c0dSLois Curfman McInnes } 814c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 815c8dd0c0dSLois Curfman McInnes } 816c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 817c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 818c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 819b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8208f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 821c8dd0c0dSLois Curfman McInnes { 822c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8236831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8246831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 825c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 826c8dd0c0dSLois Curfman McInnes } 827c8dd0c0dSLois Curfman McInnes EXTERN_C_END 828c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 82904d965bbSLois Curfman McInnes /* 8306831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 83182bf6240SBarry Smith 83204d965bbSLois Curfman McInnes Input Parameter: 83304d965bbSLois Curfman McInnes . snes - the SNES context 83404d965bbSLois Curfman McInnes 83504d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 83604d965bbSLois Curfman McInnes */ 8375615d1e5SSatish Balay #undef __FUNC__ 838b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 839f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8405e42470aSBarry Smith { 8416831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 842d132466eSBarry Smith int ierr; 8436daaf66cSBarry Smith 8443a40ed3dSBarry Smith PetscFunctionBegin; 8456831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 846d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 847d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 848d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 849d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8503a40ed3dSBarry Smith PetscFunctionReturn(0); 8515e42470aSBarry Smith } 85204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 85304d965bbSLois Curfman McInnes /* 8546831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 85504d965bbSLois Curfman McInnes 85604d965bbSLois Curfman McInnes Input Parameters: 85704d965bbSLois Curfman McInnes . SNES - the SNES context 85804d965bbSLois Curfman McInnes . viewer - visualization context 85904d965bbSLois Curfman McInnes 86004d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 86104d965bbSLois Curfman McInnes */ 8625615d1e5SSatish Balay #undef __FUNC__ 863b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 864e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 865a935fc98SLois Curfman McInnes { 8666831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 86719bcc07fSBarry Smith char *cstr; 86851695f54SLois Curfman McInnes int ierr; 8696831982aSBarry Smith PetscTruth isascii; 870a935fc98SLois Curfman McInnes 8713a40ed3dSBarry Smith PetscFunctionBegin; 8726831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8730f5bd95cSBarry Smith if (isascii) { 87419bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 87519bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 87619bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 87719bcc07fSBarry Smith else cstr = "unknown"; 878d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 879d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 8805cd90555SBarry Smith } else { 8810f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 88219bcc07fSBarry Smith } 8833a40ed3dSBarry Smith PetscFunctionReturn(0); 884a935fc98SLois Curfman McInnes } 88504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 88604d965bbSLois Curfman McInnes /* 8876831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 88804d965bbSLois Curfman McInnes 88904d965bbSLois Curfman McInnes Input Parameter: 89004d965bbSLois Curfman McInnes . snes - the SNES context 89104d965bbSLois Curfman McInnes 89204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 89304d965bbSLois Curfman McInnes */ 8945615d1e5SSatish Balay #undef __FUNC__ 895b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 896f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8975e42470aSBarry Smith { 8986831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 8995e42470aSBarry Smith char ver[16]; 9005e42470aSBarry Smith double tmp; 901f1af5d2fSBarry Smith int ierr; 902f1af5d2fSBarry Smith PetscTruth flg; 9035e42470aSBarry Smith 9043a40ed3dSBarry Smith PetscFunctionBegin; 90509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 90625cdf11fSBarry Smith if (flg) { 9075e42470aSBarry Smith ls->alpha = tmp; 9085e42470aSBarry Smith } 90909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 91025cdf11fSBarry Smith if (flg) { 9115e42470aSBarry Smith ls->maxstep = tmp; 9125e42470aSBarry Smith } 91309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 91425cdf11fSBarry Smith if (flg) { 9155e42470aSBarry Smith ls->steptol = tmp; 9165e42470aSBarry Smith } 91709d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 91825cdf11fSBarry Smith if (flg) { 919c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9200f5bd95cSBarry Smith 921c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 922c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 923c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 924c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9250f5bd95cSBarry Smith 9260f5bd95cSBarry Smith if (isbasic) { 927af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9280f5bd95cSBarry Smith } else if (isnonorms) { 929af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9300f5bd95cSBarry Smith } else if (isquad) { 931af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9320f5bd95cSBarry Smith } else if (iscubic) { 933af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9345e42470aSBarry Smith } 935a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9365e42470aSBarry Smith } 9373a40ed3dSBarry Smith PetscFunctionReturn(0); 9385e42470aSBarry Smith } 93904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 94004d965bbSLois Curfman McInnes /* 9416831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9426831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 94304d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 94404d965bbSLois Curfman McInnes 94504d965bbSLois Curfman McInnes 94604d965bbSLois Curfman McInnes Input Parameter: 94704d965bbSLois Curfman McInnes . snes - the SNES context 94804d965bbSLois Curfman McInnes 94904d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 95004d965bbSLois Curfman McInnes */ 951fb2e594dSBarry Smith EXTERN_C_BEGIN 9525615d1e5SSatish Balay #undef __FUNC__ 953b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 954f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9555e42470aSBarry Smith { 95682bf6240SBarry Smith int ierr; 9576831982aSBarry Smith SNES_EQ_LS *neP; 9585e42470aSBarry Smith 9593a40ed3dSBarry Smith PetscFunctionBegin; 960a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 961a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 962a8c6a408SBarry Smith } 96382bf6240SBarry Smith 964f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 965f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 966f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 967f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 968f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 969f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 970f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9715baf8537SBarry Smith snes->nwork = 0; 9725e42470aSBarry Smith 9736831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 9746831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 9755e42470aSBarry Smith snes->data = (void*)neP; 9765e42470aSBarry Smith neP->alpha = 1.e-4; 9775e42470aSBarry Smith neP->maxstep = 1.e8; 9785e42470aSBarry Smith neP->steptol = 1.e-12; 9795e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 980c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 981c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 982c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 98382bf6240SBarry Smith 984f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 985*0c97f09dSSatish Balay SNESSetLineSearch_LS);CHKERRQ(ierr); 986f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 987*0c97f09dSSatish Balay SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 98882bf6240SBarry Smith 9893a40ed3dSBarry Smith PetscFunctionReturn(0); 9905e42470aSBarry Smith } 991fb2e594dSBarry Smith EXTERN_C_END 99282bf6240SBarry Smith 99382bf6240SBarry Smith 99482bf6240SBarry Smith 995