1*329f5518SBarry Smith /*$Id: ls.c,v 1.149 1999/12/20 22:57:33 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__ 635615d1e5SSatish Balay #define __FUNC__ "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; 69*329f5518SBarry 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 825334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 83cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 840f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 8584c569c9SBarry Smith snes->iter = 0; 865e42470aSBarry Smith snes->norm = fnorm; 870f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 8884c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 8994a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 905e76c431SBarry Smith 91184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 9294a9d846SBarry Smith 93d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 94d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 95d034289bSBarry Smith 965e76c431SBarry Smith for (i=0; i<maxits; i++) { 975e76c431SBarry Smith 98ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 995334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1005334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 10178b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 1023505fcc1SBarry Smith snes->linear_its += lits; 103af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 104ea4d3ed3SLois Curfman McInnes 105ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 106ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 107ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 108ea4d3ed3SLois Curfman McInnes */ 10981b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 110af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 111af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 1123505fcc1SBarry Smith if (lsfail) { 1133505fcc1SBarry Smith snes->nfailures++; 1143505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 1153505fcc1SBarry Smith break; 1163505fcc1SBarry Smith } 1175e76c431SBarry Smith 11839e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 11939e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1205e76c431SBarry Smith fnorm = gnorm; 1215e76c431SBarry Smith 1220f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 12384c569c9SBarry Smith snes->iter = i+1; 1245e42470aSBarry Smith snes->norm = fnorm; 1250f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 12684c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 12794a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1285e76c431SBarry Smith 1295e76c431SBarry Smith /* Test for convergence */ 13029e0b56fSBarry Smith if (snes->converged) { 13129e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1323505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1333505fcc1SBarry Smith if (snes->reason) { 13416c95cacSBarry Smith break; 13516c95cacSBarry Smith } 13616c95cacSBarry Smith } 13729e0b56fSBarry Smith } 13839e2f89bSBarry Smith if (X != snes->vec_sol) { 139393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 14039e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 14139e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 14239e2f89bSBarry Smith } 14352392280SLois Curfman McInnes if (i == maxits) { 144981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 14552392280SLois Curfman McInnes i--; 1463505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 14752392280SLois Curfman McInnes } 1480f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1490f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1505e42470aSBarry Smith *outits = i+1; 1513a40ed3dSBarry Smith PetscFunctionReturn(0); 1525e76c431SBarry Smith } 15304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 15404d965bbSLois Curfman McInnes /* 15504d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 1566831982aSBarry Smith of the SNESEQLS nonlinear solver. 15704d965bbSLois Curfman McInnes 15804d965bbSLois Curfman McInnes Input Parameter: 15904d965bbSLois Curfman McInnes . snes - the SNES context 16004d965bbSLois Curfman McInnes . x - the solution vector 16104d965bbSLois Curfman McInnes 16204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 16304d965bbSLois Curfman McInnes 16404d965bbSLois Curfman McInnes Notes: 16504d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 16604d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 16704d965bbSLois Curfman McInnes the call to SNESSolve(). 16804d965bbSLois Curfman McInnes */ 1695615d1e5SSatish Balay #undef __FUNC__ 1705615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS" 171f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1725e76c431SBarry Smith { 1735e42470aSBarry Smith int ierr; 1743a40ed3dSBarry Smith 1753a40ed3dSBarry Smith PetscFunctionBegin; 17681b6cf68SLois Curfman McInnes snes->nwork = 4; 177d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1785e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 17981b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1803a40ed3dSBarry Smith PetscFunctionReturn(0); 1815e76c431SBarry Smith } 18204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 18304d965bbSLois Curfman McInnes /* 1846831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 18504d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 18604d965bbSLois Curfman McInnes 18704d965bbSLois Curfman McInnes Input Parameter: 18804d965bbSLois Curfman McInnes . snes - the SNES context 18904d965bbSLois Curfman McInnes 190de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 19104d965bbSLois Curfman McInnes */ 1925615d1e5SSatish Balay #undef __FUNC__ 193d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS" 194e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1955e76c431SBarry Smith { 196393d2d9aSLois Curfman McInnes int ierr; 1973a40ed3dSBarry Smith 1983a40ed3dSBarry Smith PetscFunctionBegin; 1995baf8537SBarry Smith if (snes->nwork) { 2004b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 20121c89e3eSBarry Smith } 2025c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 2045e76c431SBarry Smith } 20504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2065615d1e5SSatish Balay #undef __FUNC__ 2075615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch" 20882bf6240SBarry Smith 2094b828684SBarry Smith /*@C 2105e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2115e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2125e76c431SBarry Smith to serve as a template and is not recommended for general use. 2135e76c431SBarry Smith 214c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 215c7afd0dbSLois Curfman McInnes 2165e76c431SBarry Smith Input Parameters: 217c7afd0dbSLois Curfman McInnes + snes - nonlinear context 218af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2195e76c431SBarry Smith . x - current iterate 2205e76c431SBarry Smith . f - residual evaluated at x 2215e76c431SBarry Smith . y - search direction (contains new iterate on output) 2225e76c431SBarry Smith . w - work vector 223c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2245e76c431SBarry Smith 225c4a48953SLois Curfman McInnes Output Parameters: 226c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2275e76c431SBarry Smith . y - new iterate (contains search direction on input) 2285e76c431SBarry Smith . gnorm - 2-norm of g 2295e76c431SBarry Smith . ynorm - 2-norm of search length 230c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 231fee21e36SBarry Smith 232c4a48953SLois Curfman McInnes Options Database Key: 233c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 234c4a48953SLois Curfman McInnes 23536851e7fSLois Curfman McInnes Level: advanced 23636851e7fSLois Curfman McInnes 23728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 23828ae5a21SLois Curfman McInnes 239f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 240af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2415e76c431SBarry Smith @*/ 242af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 243*329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 2445e76c431SBarry Smith { 2455e42470aSBarry Smith int ierr; 2465334005bSBarry Smith Scalar mone = -1.0; 2476831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 2488f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2495334005bSBarry Smith 2503a40ed3dSBarry Smith PetscFunctionBegin; 251761aaf1bSLois Curfman McInnes *flag = 0; 2527857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 253a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 254ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 2558f99978dSLois Curfman McInnes if (neP->CheckStep) { 2568f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 2578f99978dSLois Curfman McInnes } 258ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 259a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 2607857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2613a40ed3dSBarry Smith PetscFunctionReturn(0); 2625e76c431SBarry Smith } 26304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2645615d1e5SSatish Balay #undef __FUNC__ 26529e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 26682bf6240SBarry Smith 26729e0b56fSBarry Smith /*@C 26829e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 26929e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 27029e0b56fSBarry Smith even compute the norm of the function or search direction; this 27129e0b56fSBarry Smith is intended only when you know the full step is fine and are 27229e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 273c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 274c7afd0dbSLois Curfman McInnes 275c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 27629e0b56fSBarry Smith 27729e0b56fSBarry Smith Input Parameters: 278c7afd0dbSLois Curfman McInnes + snes - nonlinear context 279af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 28029e0b56fSBarry Smith . x - current iterate 28129e0b56fSBarry Smith . f - residual evaluated at x 28229e0b56fSBarry Smith . y - search direction (contains new iterate on output) 28329e0b56fSBarry Smith . w - work vector 284c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 28529e0b56fSBarry Smith 28629e0b56fSBarry Smith Output Parameters: 287c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 28829e0b56fSBarry Smith . gnorm - not changed 28929e0b56fSBarry Smith . ynorm - not changed 290c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 291fee21e36SBarry Smith 29229e0b56fSBarry Smith Options Database Key: 293c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 29429e0b56fSBarry Smith 2958cbba510SBarry Smith Notes: 296ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 297ea56f5baSLois Curfman McInnes either the options 298ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 299ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 300*329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 301*329f5518SBarry Smith otherwise, the SNES solver will generate an error. 302*329f5518SBarry Smith 303*329f5518SBarry Smith During the final iteration this will not evaluate the function at 304*329f5518SBarry Smith the solution point. This is to save a function evaluation while 305*329f5518SBarry Smith using pseudo-timestepping. 3068cbba510SBarry Smith 307ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 308ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 309ea56f5baSLois Curfman McInnes correct, since they are not computed. 310ea56f5baSLois Curfman McInnes 311ea56f5baSLois Curfman McInnes Level: advanced 3128cbba510SBarry Smith 31329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 31429e0b56fSBarry Smith 31529e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 31629e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 31729e0b56fSBarry Smith @*/ 318af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 319*329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 32029e0b56fSBarry Smith { 32129e0b56fSBarry Smith int ierr; 32229e0b56fSBarry Smith Scalar mone = -1.0; 3236831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3248f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 32529e0b56fSBarry Smith 3263a40ed3dSBarry Smith PetscFunctionBegin; 3278cbba510SBarry Smith *flag = 0; 32829e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 32929e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3308f99978dSLois Curfman McInnes if (neP->CheckStep) { 3318f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3328f99978dSLois Curfman McInnes } 333*329f5518SBarry Smith 334*329f5518SBarry Smith /* don't evaluate function the last time through */ 335*329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 33629e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 337*329f5518SBarry Smith } 33829e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 34029e0b56fSBarry Smith } 34104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34229e0b56fSBarry Smith #undef __FUNC__ 3435615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3444b828684SBarry Smith /*@C 345f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3465e76c431SBarry Smith 347c7afd0dbSLois Curfman McInnes Collective on SNES 348c7afd0dbSLois Curfman McInnes 3495e76c431SBarry Smith Input Parameters: 350c7afd0dbSLois Curfman McInnes + snes - nonlinear context 351af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3525e76c431SBarry Smith . x - current iterate 3535e76c431SBarry Smith . f - residual evaluated at x 3545e76c431SBarry Smith . y - search direction (contains new iterate on output) 3555e76c431SBarry Smith . w - work vector 356c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3575e76c431SBarry Smith 358393d2d9aSLois Curfman McInnes Output Parameters: 359c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3605e76c431SBarry Smith . y - new iterate (contains search direction on input) 3615e76c431SBarry Smith . gnorm - 2-norm of g 3625e76c431SBarry Smith . ynorm - 2-norm of search length 363c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 364fee21e36SBarry Smith 365c4a48953SLois Curfman McInnes Options Database Key: 366c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 367c4a48953SLois Curfman McInnes 3685e76c431SBarry Smith Notes: 3695e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3705e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3715e76c431SBarry Smith 37236851e7fSLois Curfman McInnes Level: advanced 37336851e7fSLois Curfman McInnes 37428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 37528ae5a21SLois Curfman McInnes 376af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3775e76c431SBarry Smith @*/ 378af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 379*329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3805e76c431SBarry Smith { 381406556e6SLois Curfman McInnes /* 382406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 383406556e6SLois Curfman McInnes minimization problem: 384406556e6SLois Curfman McInnes min z(x): R^n -> R, 385406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 386406556e6SLois Curfman McInnes */ 387406556e6SLois Curfman McInnes 388*329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 389*329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 390aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 3915e42470aSBarry Smith Scalar cinitslope,clambda; 3926b5873e3SBarry Smith #endif 3935e42470aSBarry Smith int ierr,count; 3946831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3955334005bSBarry Smith Scalar mone = -1.0,scale; 3968f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3975e76c431SBarry Smith 3983a40ed3dSBarry Smith PetscFunctionBegin; 3997857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 400761aaf1bSLois Curfman McInnes *flag = 0; 4015e76c431SBarry Smith alpha = neP->alpha; 4025e76c431SBarry Smith maxstep = neP->maxstep; 4035e76c431SBarry Smith steptol = neP->steptol; 4045e76c431SBarry Smith 405cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 406a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 407a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 408a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 409a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 410a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 411ad922ac9SBarry Smith goto theend1; 41294a9d846SBarry Smith } 4135e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4145e42470aSBarry Smith scale = maxstep/(*ynorm); 415aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 416*329f5518SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscRealPart(scale)); 4176b5873e3SBarry Smith #else 41894a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 4196b5873e3SBarry Smith #endif 420393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4215e76c431SBarry Smith *ynorm = maxstep; 4225e76c431SBarry Smith } 4235e76c431SBarry Smith minlambda = steptol/(*ynorm); 424a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 425aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 426a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 427*329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 4285e42470aSBarry Smith #else 429a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4305e42470aSBarry Smith #endif 4315e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4325e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4335e76c431SBarry Smith 434393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4355334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 43678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 437313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 438406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 439393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 44094a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 441ad922ac9SBarry Smith goto theend1; 4425e76c431SBarry Smith } 4435e76c431SBarry Smith 4445e76c431SBarry Smith /* Fit points with quadratic */ 445313b5538SBarry Smith lambda = 1.0; 446a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4475e76c431SBarry Smith lambdaprev = lambda; 4485e76c431SBarry Smith gnormprev = *gnorm; 449*329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 450ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 451ddd12767SBarry Smith else lambda = lambdatemp; 452393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 453ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 454aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 455ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4565e42470aSBarry Smith #else 457ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4585e42470aSBarry Smith #endif 45978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 460cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 461406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 462393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 463ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 464ad922ac9SBarry Smith goto theend1; 4655e76c431SBarry Smith } 4665e76c431SBarry Smith 4675e76c431SBarry Smith /* Fit points with cubic */ 4685e76c431SBarry Smith count = 1; 4695e76c431SBarry Smith while (1) { 4705e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4712b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4722b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 473ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 474393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 475761aaf1bSLois Curfman McInnes *flag = -1; break; 4765e76c431SBarry Smith } 477406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 478406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 479ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4802b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4815e76c431SBarry Smith d = b*b - 3*a*initslope; 4825e76c431SBarry Smith if (d < 0.0) d = 0.0; 4835e76c431SBarry Smith if (a == 0.0) { 4845e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4855e76c431SBarry Smith } else { 4865e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4875e76c431SBarry Smith } 4885e76c431SBarry Smith lambdaprev = lambda; 4895e76c431SBarry Smith gnormprev = *gnorm; 490*329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 491*329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 4925e76c431SBarry Smith else lambda = lambdatemp; 493393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 494ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 495aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 496ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 497393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4985e42470aSBarry Smith #else 499ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5005e42470aSBarry Smith #endif 50178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 502cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 503406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 504393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 505ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5062715565aSLois Curfman McInnes goto theend1; 5072b022350SLois Curfman McInnes } else { 5082b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5095e76c431SBarry Smith } 5105e76c431SBarry Smith count++; 5115e76c431SBarry Smith } 512ad922ac9SBarry Smith theend1: 5138f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5148f99978dSLois Curfman McInnes if (neP->CheckStep) { 5158f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5168f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5178f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5188f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5198f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5208f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5218f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5228f99978dSLois Curfman McInnes } 5238f99978dSLois Curfman McInnes } 5247857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5253a40ed3dSBarry Smith PetscFunctionReturn(0); 5265e76c431SBarry Smith } 52704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5285615d1e5SSatish Balay #undef __FUNC__ 5295615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 5304b828684SBarry Smith /*@C 531f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5325e76c431SBarry Smith 533c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 534c7afd0dbSLois Curfman McInnes 5355e42470aSBarry Smith Input Parameters: 536c7afd0dbSLois Curfman McInnes + snes - the SNES context 537af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5385e42470aSBarry Smith . x - current iterate 5395e42470aSBarry Smith . f - residual evaluated at x 5405e42470aSBarry Smith . y - search direction (contains new iterate on output) 5415e42470aSBarry Smith . w - work vector 542c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5435e42470aSBarry Smith 544c4a48953SLois Curfman McInnes Output Parameters: 545c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5465e42470aSBarry Smith . y - new iterate (contains search direction on input) 5475e42470aSBarry Smith . gnorm - 2-norm of g 5485e42470aSBarry Smith . ynorm - 2-norm of search length 549c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 550fee21e36SBarry Smith 551c4a48953SLois Curfman McInnes Options Database Key: 552c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 553c4a48953SLois Curfman McInnes 5545e42470aSBarry Smith Notes: 5556831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 5565e42470aSBarry Smith 55736851e7fSLois Curfman McInnes Level: advanced 55836851e7fSLois Curfman McInnes 559f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 560f59ffdedSLois Curfman McInnes 561af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5625e42470aSBarry Smith @*/ 563af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 564*329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 5655e76c431SBarry Smith { 566406556e6SLois Curfman McInnes /* 567406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 568406556e6SLois Curfman McInnes minimization problem: 569406556e6SLois Curfman McInnes min z(x): R^n -> R, 570406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 571406556e6SLois Curfman McInnes */ 572*329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 573aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5745e42470aSBarry Smith Scalar cinitslope,clambda; 5756b5873e3SBarry Smith #endif 5765e42470aSBarry Smith int ierr,count; 5776831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 5785334005bSBarry Smith Scalar mone = -1.0,scale; 5798f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5805e76c431SBarry Smith 5813a40ed3dSBarry Smith PetscFunctionBegin; 5827857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 583761aaf1bSLois Curfman McInnes *flag = 0; 5845e42470aSBarry Smith alpha = neP->alpha; 5855e42470aSBarry Smith maxstep = neP->maxstep; 5865e42470aSBarry Smith steptol = neP->steptol; 5875e76c431SBarry Smith 5883505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 589b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 59094a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 591b37302c6SLois Curfman McInnes *gnorm = fnorm; 592b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 593b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 594ad922ac9SBarry Smith goto theend2; 59594a9d846SBarry Smith } 5965e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5975e42470aSBarry Smith scale = maxstep/(*ynorm); 598393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5995e42470aSBarry Smith *ynorm = maxstep; 6005e76c431SBarry Smith } 6015e42470aSBarry Smith minlambda = steptol/(*ynorm); 602a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 603aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 604a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 605*329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6065e42470aSBarry Smith #else 607a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6085e42470aSBarry Smith #endif 6095e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6105e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6115e42470aSBarry Smith 612393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6135334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 61478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 615cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 616406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 617393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 61894a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 619ad922ac9SBarry Smith goto theend2; 6205e42470aSBarry Smith } 6215e42470aSBarry Smith 6225e42470aSBarry Smith /* Fit points with quadratic */ 623313b5538SBarry Smith lambda = 1.0; 6245e42470aSBarry Smith count = 1; 6255e42470aSBarry Smith while (1) { 6265e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 627981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 628981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 629ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 630393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 631761aaf1bSLois Curfman McInnes *flag = -1; break; 6325e42470aSBarry Smith } 633a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 634*329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 635*329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 636*329f5518SBarry Smith else lambda = lambdatemp; 637393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6383505fcc1SBarry Smith lambdaneg = -lambda; 639aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6403505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6415e42470aSBarry Smith #else 6423505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6435e42470aSBarry Smith #endif 64478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 645cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 646406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 647393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 648981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 64906259719SBarry Smith break; 6505e42470aSBarry Smith } 6515e42470aSBarry Smith count++; 6525e42470aSBarry Smith } 653ad922ac9SBarry Smith theend2: 6548f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6558f99978dSLois Curfman McInnes if (neP->CheckStep) { 6568f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6578f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6588f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6598f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6608f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6618f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6628f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6638f99978dSLois Curfman McInnes } 6648f99978dSLois Curfman McInnes } 6657857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6663a40ed3dSBarry Smith PetscFunctionReturn(0); 6675e76c431SBarry Smith } 66804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6695615d1e5SSatish Balay #undef __FUNC__ 670d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 671c9e6a524SLois Curfman McInnes /*@C 67277c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 6736831982aSBarry Smith by the method SNESEQLS. 6745e76c431SBarry Smith 6755e76c431SBarry Smith Input Parameters: 676c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 677af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 678c7afd0dbSLois Curfman McInnes - func - pointer to int function 6795e76c431SBarry Smith 680fee21e36SBarry Smith Collective on SNES 681fee21e36SBarry Smith 682c4a48953SLois Curfman McInnes Available Routines: 683c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 684f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 685f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 686af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6875e76c431SBarry Smith 688c4a48953SLois Curfman McInnes Options Database Keys: 689af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 690c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 691c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 692c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 693c4a48953SLois Curfman McInnes 6945e76c431SBarry Smith Calling sequence of func: 695bd208895SLois Curfman McInnes .vb 696af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 697*329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 698bd208895SLois Curfman McInnes .ve 6995e76c431SBarry Smith 7005e76c431SBarry Smith Input parameters for func: 701c7afd0dbSLois Curfman McInnes + snes - nonlinear context 702af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7035e76c431SBarry Smith . x - current iterate 7045e76c431SBarry Smith . f - residual evaluated at x 7055e76c431SBarry Smith . y - search direction (contains new iterate on output) 7065e76c431SBarry Smith . w - work vector 707c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7085e76c431SBarry Smith 7095e76c431SBarry Smith Output parameters for func: 710c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7115e76c431SBarry Smith . y - new iterate (contains search direction on input) 7125e76c431SBarry Smith . gnorm - 2-norm of g 7135e76c431SBarry Smith . ynorm - 2-norm of search length 714c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 715761aaf1bSLois Curfman McInnes on failure. 716f59ffdedSLois Curfman McInnes 71736851e7fSLois Curfman McInnes Level: advanced 71836851e7fSLois Curfman McInnes 719f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 720f59ffdedSLois Curfman McInnes 7214ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7225e76c431SBarry Smith @*/ 723*329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7245e76c431SBarry Smith { 725*329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 72682bf6240SBarry Smith 7273a40ed3dSBarry Smith PetscFunctionBegin; 72894d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 72982bf6240SBarry Smith if (f) { 730af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 73182bf6240SBarry Smith } 7323a40ed3dSBarry Smith PetscFunctionReturn(0); 7335e76c431SBarry Smith } 73404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 735fb2e594dSBarry Smith EXTERN_C_BEGIN 73682bf6240SBarry Smith #undef __FUNC__ 73782bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 738af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 739af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 74082bf6240SBarry Smith { 74182bf6240SBarry Smith PetscFunctionBegin; 7426831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7436831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 74482bf6240SBarry Smith PetscFunctionReturn(0); 74582bf6240SBarry Smith } 746fb2e594dSBarry Smith EXTERN_C_END 74704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 748c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 749c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck" 750c8dd0c0dSLois Curfman McInnes /*@C 751530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7526831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 753c8dd0c0dSLois Curfman McInnes 754c8dd0c0dSLois Curfman McInnes Input Parameters: 755c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 756c8dd0c0dSLois Curfman McInnes . func - pointer to int function 757c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 758c8dd0c0dSLois Curfman McInnes 759c8dd0c0dSLois Curfman McInnes Collective on SNES 760c8dd0c0dSLois Curfman McInnes 761c8dd0c0dSLois Curfman McInnes Calling sequence of func: 762c8dd0c0dSLois Curfman McInnes .vb 763b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 764c8dd0c0dSLois Curfman McInnes .ve 765b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 766b1ae27eaSLois Curfman McInnes on failure. 767c8dd0c0dSLois Curfman McInnes 768c8dd0c0dSLois Curfman McInnes Input parameters for func: 769c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 770c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7718f99978dSLois Curfman McInnes - x - current candidate iterate 772c8dd0c0dSLois Curfman McInnes 773c8dd0c0dSLois Curfman McInnes Output parameters for func: 774c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 775c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 776c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 777c8dd0c0dSLois Curfman McInnes 778c8dd0c0dSLois Curfman McInnes Level: advanced 779c8dd0c0dSLois Curfman McInnes 7808f99978dSLois Curfman McInnes Notes: 781b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 782b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 783b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 784b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 785ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 7868f99978dSLois Curfman McInnes 787b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 788b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 789b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 790ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 791ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 792ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 793ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 794b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 7958f99978dSLois Curfman McInnes 796c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 797c8dd0c0dSLois Curfman McInnes 798c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 799c8dd0c0dSLois Curfman McInnes @*/ 8008f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 801c8dd0c0dSLois Curfman McInnes { 8028f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 803c8dd0c0dSLois Curfman McInnes 804c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 805c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 806c8dd0c0dSLois Curfman McInnes if (f) { 807c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 808c8dd0c0dSLois Curfman McInnes } 809c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 810c8dd0c0dSLois Curfman McInnes } 811c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 812c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 813c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 814c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS" 8158f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 816c8dd0c0dSLois Curfman McInnes { 817c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8186831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8196831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 820c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 821c8dd0c0dSLois Curfman McInnes } 822c8dd0c0dSLois Curfman McInnes EXTERN_C_END 823c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 82404d965bbSLois Curfman McInnes /* 8256831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 82682bf6240SBarry Smith 82704d965bbSLois Curfman McInnes Input Parameter: 82804d965bbSLois Curfman McInnes . snes - the SNES context 82904d965bbSLois Curfman McInnes 83004d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 83104d965bbSLois Curfman McInnes */ 8325615d1e5SSatish Balay #undef __FUNC__ 833d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 834f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8355e42470aSBarry Smith { 8366831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 837d132466eSBarry Smith int ierr; 8386daaf66cSBarry Smith 8393a40ed3dSBarry Smith PetscFunctionBegin; 8406831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 841d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 842d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 843d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 844d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8453a40ed3dSBarry Smith PetscFunctionReturn(0); 8465e42470aSBarry Smith } 84704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 84804d965bbSLois Curfman McInnes /* 8496831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 85004d965bbSLois Curfman McInnes 85104d965bbSLois Curfman McInnes Input Parameters: 85204d965bbSLois Curfman McInnes . SNES - the SNES context 85304d965bbSLois Curfman McInnes . viewer - visualization context 85404d965bbSLois Curfman McInnes 85504d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 85604d965bbSLois Curfman McInnes */ 8575615d1e5SSatish Balay #undef __FUNC__ 858d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 859e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 860a935fc98SLois Curfman McInnes { 8616831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 86219bcc07fSBarry Smith char *cstr; 86351695f54SLois Curfman McInnes int ierr; 8646831982aSBarry Smith PetscTruth isascii; 865a935fc98SLois Curfman McInnes 8663a40ed3dSBarry Smith PetscFunctionBegin; 8676831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8680f5bd95cSBarry Smith if (isascii) { 86919bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 87019bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 87119bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 87219bcc07fSBarry Smith else cstr = "unknown"; 873d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 874d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 8755cd90555SBarry Smith } else { 8760f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 87719bcc07fSBarry Smith } 8783a40ed3dSBarry Smith PetscFunctionReturn(0); 879a935fc98SLois Curfman McInnes } 88004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 88104d965bbSLois Curfman McInnes /* 8826831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 88304d965bbSLois Curfman McInnes 88404d965bbSLois Curfman McInnes Input Parameter: 88504d965bbSLois Curfman McInnes . snes - the SNES context 88604d965bbSLois Curfman McInnes 88704d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 88804d965bbSLois Curfman McInnes */ 8895615d1e5SSatish Balay #undef __FUNC__ 8905615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 891f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8925e42470aSBarry Smith { 8936831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 8945e42470aSBarry Smith char ver[16]; 8955e42470aSBarry Smith double tmp; 896f1af5d2fSBarry Smith int ierr; 897f1af5d2fSBarry Smith PetscTruth flg; 8985e42470aSBarry Smith 8993a40ed3dSBarry Smith PetscFunctionBegin; 90009d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 90125cdf11fSBarry Smith if (flg) { 9025e42470aSBarry Smith ls->alpha = tmp; 9035e42470aSBarry Smith } 90409d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 90525cdf11fSBarry Smith if (flg) { 9065e42470aSBarry Smith ls->maxstep = tmp; 9075e42470aSBarry Smith } 90809d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 90925cdf11fSBarry Smith if (flg) { 9105e42470aSBarry Smith ls->steptol = tmp; 9115e42470aSBarry Smith } 91209d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 91325cdf11fSBarry Smith if (flg) { 914c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9150f5bd95cSBarry Smith 916c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 917c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 918c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 919c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9200f5bd95cSBarry Smith 9210f5bd95cSBarry Smith if (isbasic) { 922af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9230f5bd95cSBarry Smith } else if (isnonorms) { 924af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9250f5bd95cSBarry Smith } else if (isquad) { 926af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9270f5bd95cSBarry Smith } else if (iscubic) { 928af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9295e42470aSBarry Smith } 930a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9315e42470aSBarry Smith } 9323a40ed3dSBarry Smith PetscFunctionReturn(0); 9335e42470aSBarry Smith } 93404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 93504d965bbSLois Curfman McInnes /* 9366831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9376831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 93804d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 93904d965bbSLois Curfman McInnes 94004d965bbSLois Curfman McInnes 94104d965bbSLois Curfman McInnes Input Parameter: 94204d965bbSLois Curfman McInnes . snes - the SNES context 94304d965bbSLois Curfman McInnes 94404d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 94504d965bbSLois Curfman McInnes */ 946fb2e594dSBarry Smith EXTERN_C_BEGIN 9475615d1e5SSatish Balay #undef __FUNC__ 9485615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 949f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9505e42470aSBarry Smith { 95182bf6240SBarry Smith int ierr; 9526831982aSBarry Smith SNES_EQ_LS *neP; 9535e42470aSBarry Smith 9543a40ed3dSBarry Smith PetscFunctionBegin; 955a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 956a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 957a8c6a408SBarry Smith } 95882bf6240SBarry Smith 959f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 960f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 961f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 962f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 963f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 964f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 965f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9665baf8537SBarry Smith snes->nwork = 0; 9675e42470aSBarry Smith 9686831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 9696831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 9705e42470aSBarry Smith snes->data = (void*)neP; 9715e42470aSBarry Smith neP->alpha = 1.e-4; 9725e42470aSBarry Smith neP->maxstep = 1.e8; 9735e42470aSBarry Smith neP->steptol = 1.e-12; 9745e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 975c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 976c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 977c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 97882bf6240SBarry Smith 979f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 980e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 981f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 982c8dd0c0dSLois Curfman McInnes (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 98382bf6240SBarry Smith 9843a40ed3dSBarry Smith PetscFunctionReturn(0); 9855e42470aSBarry Smith } 986fb2e594dSBarry Smith EXTERN_C_END 98782bf6240SBarry Smith 98882bf6240SBarry Smith 98982bf6240SBarry Smith 990