1*c38d4ed2SBarry Smith /*$Id: ls.c,v 1.147 1999/11/05 14:47:10 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; 6984c569c9SBarry Smith double fnorm, gnorm, xnorm, ynorm; 705e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 71184914b5SBarry Smith SNESConvergedReason reason; 725e76c431SBarry Smith 733a40ed3dSBarry Smith PetscFunctionBegin; 74184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 75184914b5SBarry Smith 765e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 775e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 7839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 795e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 805e42470aSBarry Smith G = snes->work[1]; 815e42470aSBarry Smith W = snes->work[2]; 825e76c431SBarry Smith 835334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 84cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 850f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 8684c569c9SBarry Smith snes->iter = 0; 875e42470aSBarry Smith snes->norm = fnorm; 880f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 8984c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 9094a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 915e76c431SBarry Smith 92184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 9394a9d846SBarry Smith 94d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 95d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 96d034289bSBarry Smith 975e76c431SBarry Smith for ( i=0; i<maxits; i++ ) { 985e76c431SBarry Smith 99ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1005334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1015334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 10278b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 1037a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(lits); 104af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 105ea4d3ed3SLois Curfman McInnes 106ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 107ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 108ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 109ea4d3ed3SLois Curfman McInnes */ 11081b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 111af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 112af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 113761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 1145e76c431SBarry Smith 11539e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 11639e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1175e76c431SBarry Smith fnorm = gnorm; 1185e76c431SBarry Smith 1190f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 12084c569c9SBarry Smith snes->iter = i+1; 1215e42470aSBarry Smith snes->norm = fnorm; 1220f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 12384c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 12494a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1255e76c431SBarry Smith 1265e76c431SBarry Smith /* Test for convergence */ 12729e0b56fSBarry Smith if (snes->converged) { 12829e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 129184914b5SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 130184914b5SBarry Smith if (reason) { 13116c95cacSBarry Smith break; 13216c95cacSBarry Smith } 13316c95cacSBarry Smith } 13429e0b56fSBarry Smith } 13539e2f89bSBarry Smith if (X != snes->vec_sol) { 136393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 13739e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 13839e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 13939e2f89bSBarry Smith } 14052392280SLois Curfman McInnes if (i == maxits) { 141981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 14252392280SLois Curfman McInnes i--; 143184914b5SBarry Smith reason = SNES_DIVERGED_MAX_IT; 14452392280SLois Curfman McInnes } 1450f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 146184914b5SBarry Smith snes->reason = reason; 1470f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1485e42470aSBarry Smith *outits = i+1; 1493a40ed3dSBarry Smith PetscFunctionReturn(0); 1505e76c431SBarry Smith } 15104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 15204d965bbSLois Curfman McInnes /* 15304d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 1546831982aSBarry Smith of the SNESEQLS nonlinear solver. 15504d965bbSLois Curfman McInnes 15604d965bbSLois Curfman McInnes Input Parameter: 15704d965bbSLois Curfman McInnes . snes - the SNES context 15804d965bbSLois Curfman McInnes . x - the solution vector 15904d965bbSLois Curfman McInnes 16004d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 16104d965bbSLois Curfman McInnes 16204d965bbSLois Curfman McInnes Notes: 16304d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 16404d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 16504d965bbSLois Curfman McInnes the call to SNESSolve(). 16604d965bbSLois Curfman McInnes */ 1675615d1e5SSatish Balay #undef __FUNC__ 1685615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS" 169f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1705e76c431SBarry Smith { 1715e42470aSBarry Smith int ierr; 1723a40ed3dSBarry Smith 1733a40ed3dSBarry Smith PetscFunctionBegin; 17481b6cf68SLois Curfman McInnes snes->nwork = 4; 175d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1765e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 17781b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1783a40ed3dSBarry Smith PetscFunctionReturn(0); 1795e76c431SBarry Smith } 18004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 18104d965bbSLois Curfman McInnes /* 1826831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 18304d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 18404d965bbSLois Curfman McInnes 18504d965bbSLois Curfman McInnes Input Parameter: 18604d965bbSLois Curfman McInnes . snes - the SNES context 18704d965bbSLois Curfman McInnes 188de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 18904d965bbSLois Curfman McInnes */ 1905615d1e5SSatish Balay #undef __FUNC__ 191d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS" 192e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1935e76c431SBarry Smith { 194393d2d9aSLois Curfman McInnes int ierr; 1953a40ed3dSBarry Smith 1963a40ed3dSBarry Smith PetscFunctionBegin; 1975baf8537SBarry Smith if (snes->nwork) { 1984b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 19921c89e3eSBarry Smith } 2005c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2013a40ed3dSBarry Smith PetscFunctionReturn(0); 2025e76c431SBarry Smith } 20304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2045615d1e5SSatish Balay #undef __FUNC__ 2055615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch" 20682bf6240SBarry Smith 2074b828684SBarry Smith /*@C 2085e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2095e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2105e76c431SBarry Smith to serve as a template and is not recommended for general use. 2115e76c431SBarry Smith 212c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 213c7afd0dbSLois Curfman McInnes 2145e76c431SBarry Smith Input Parameters: 215c7afd0dbSLois Curfman McInnes + snes - nonlinear context 216af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2175e76c431SBarry Smith . x - current iterate 2185e76c431SBarry Smith . f - residual evaluated at x 2195e76c431SBarry Smith . y - search direction (contains new iterate on output) 2205e76c431SBarry Smith . w - work vector 221c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2225e76c431SBarry Smith 223c4a48953SLois Curfman McInnes Output Parameters: 224c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2255e76c431SBarry Smith . y - new iterate (contains search direction on input) 2265e76c431SBarry Smith . gnorm - 2-norm of g 2275e76c431SBarry Smith . ynorm - 2-norm of search length 228c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 229fee21e36SBarry Smith 230c4a48953SLois Curfman McInnes Options Database Key: 231c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 232c4a48953SLois Curfman McInnes 23336851e7fSLois Curfman McInnes Level: advanced 23436851e7fSLois Curfman McInnes 23528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 23628ae5a21SLois Curfman McInnes 237f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 238af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2395e76c431SBarry Smith @*/ 240af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 241761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 2425e76c431SBarry Smith { 2435e42470aSBarry Smith int ierr; 2445334005bSBarry Smith Scalar mone = -1.0; 2456831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data; 2468f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2475334005bSBarry Smith 2483a40ed3dSBarry Smith PetscFunctionBegin; 249761aaf1bSLois Curfman McInnes *flag = 0; 2507857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 251a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 252ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 2538f99978dSLois Curfman McInnes if (neP->CheckStep) { 2548f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 2558f99978dSLois Curfman McInnes } 256ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 257a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 2587857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2593a40ed3dSBarry Smith PetscFunctionReturn(0); 2605e76c431SBarry Smith } 26104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2625615d1e5SSatish Balay #undef __FUNC__ 26329e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 26482bf6240SBarry Smith 26529e0b56fSBarry Smith /*@C 26629e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 26729e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 26829e0b56fSBarry Smith even compute the norm of the function or search direction; this 26929e0b56fSBarry Smith is intended only when you know the full step is fine and are 27029e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 271c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 272c7afd0dbSLois Curfman McInnes 273c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 27429e0b56fSBarry Smith 27529e0b56fSBarry Smith Input Parameters: 276c7afd0dbSLois Curfman McInnes + snes - nonlinear context 277af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 27829e0b56fSBarry Smith . x - current iterate 27929e0b56fSBarry Smith . f - residual evaluated at x 28029e0b56fSBarry Smith . y - search direction (contains new iterate on output) 28129e0b56fSBarry Smith . w - work vector 282c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 28329e0b56fSBarry Smith 28429e0b56fSBarry Smith Output Parameters: 285c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 28629e0b56fSBarry Smith . gnorm - not changed 28729e0b56fSBarry Smith . ynorm - not changed 288c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 289fee21e36SBarry Smith 29029e0b56fSBarry Smith Options Database Key: 291c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 29229e0b56fSBarry Smith 2938cbba510SBarry Smith Notes: 294ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 295ea56f5baSLois Curfman McInnes either the options 296ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 297ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 298ea56f5baSLois Curfman McInnes SNESSetConvergenceTest(); otherwise, the SNES solver will 299ea56f5baSLois Curfman McInnes generate an error. 3008cbba510SBarry Smith 301ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 302ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 303ea56f5baSLois Curfman McInnes correct, since they are not computed. 304ea56f5baSLois Curfman McInnes 305ea56f5baSLois Curfman McInnes Level: advanced 3068cbba510SBarry Smith 30729e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 30829e0b56fSBarry Smith 30929e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 31029e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 31129e0b56fSBarry Smith @*/ 312af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 31329e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 31429e0b56fSBarry Smith { 31529e0b56fSBarry Smith int ierr; 31629e0b56fSBarry Smith Scalar mone = -1.0; 3176831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data; 3188f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 31929e0b56fSBarry Smith 3203a40ed3dSBarry Smith PetscFunctionBegin; 3218cbba510SBarry Smith *flag = 0; 32229e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 32329e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3248f99978dSLois Curfman McInnes if (neP->CheckStep) { 3258f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3268f99978dSLois Curfman McInnes } 32729e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 32829e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3293a40ed3dSBarry Smith PetscFunctionReturn(0); 33029e0b56fSBarry Smith } 33104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 33229e0b56fSBarry Smith #undef __FUNC__ 3335615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3344b828684SBarry Smith /*@C 335f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3365e76c431SBarry Smith 337c7afd0dbSLois Curfman McInnes Collective on SNES 338c7afd0dbSLois Curfman McInnes 3395e76c431SBarry Smith Input Parameters: 340c7afd0dbSLois Curfman McInnes + snes - nonlinear context 341af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3425e76c431SBarry Smith . x - current iterate 3435e76c431SBarry Smith . f - residual evaluated at x 3445e76c431SBarry Smith . y - search direction (contains new iterate on output) 3455e76c431SBarry Smith . w - work vector 346c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3475e76c431SBarry Smith 348393d2d9aSLois Curfman McInnes Output Parameters: 349c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3505e76c431SBarry Smith . y - new iterate (contains search direction on input) 3515e76c431SBarry Smith . gnorm - 2-norm of g 3525e76c431SBarry Smith . ynorm - 2-norm of search length 353c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 354fee21e36SBarry Smith 355c4a48953SLois Curfman McInnes Options Database Key: 356c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 357c4a48953SLois Curfman McInnes 3585e76c431SBarry Smith Notes: 3595e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3605e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3615e76c431SBarry Smith 36236851e7fSLois Curfman McInnes Level: advanced 36336851e7fSLois Curfman McInnes 36428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 36528ae5a21SLois Curfman McInnes 366af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3675e76c431SBarry Smith @*/ 368af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 369761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3705e76c431SBarry Smith { 371406556e6SLois Curfman McInnes /* 372406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 373406556e6SLois Curfman McInnes minimization problem: 374406556e6SLois Curfman McInnes min z(x): R^n -> R, 375406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 376406556e6SLois Curfman McInnes */ 377406556e6SLois Curfman McInnes 378ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 379ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 380aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 3815e42470aSBarry Smith Scalar cinitslope, clambda; 3826b5873e3SBarry Smith #endif 3835e42470aSBarry Smith int ierr, count; 3846831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data; 3855334005bSBarry Smith Scalar mone = -1.0, scale; 3868f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3875e76c431SBarry Smith 3883a40ed3dSBarry Smith PetscFunctionBegin; 3897857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 390761aaf1bSLois Curfman McInnes *flag = 0; 3915e76c431SBarry Smith alpha = neP->alpha; 3925e76c431SBarry Smith maxstep = neP->maxstep; 3935e76c431SBarry Smith steptol = neP->steptol; 3945e76c431SBarry Smith 395cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 396a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 397a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 398a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 399a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 400a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 401ad922ac9SBarry Smith goto theend1; 40294a9d846SBarry Smith } 4035e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4045e42470aSBarry Smith scale = maxstep/(*ynorm); 405aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 40692318cfeSSatish Balay PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale)); 4076b5873e3SBarry Smith #else 40894a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 4096b5873e3SBarry Smith #endif 410393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4115e76c431SBarry Smith *ynorm = maxstep; 4125e76c431SBarry Smith } 4135e76c431SBarry Smith minlambda = steptol/(*ynorm); 414a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 415aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 416a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 41792318cfeSSatish Balay initslope = PetscReal(cinitslope); 4185e42470aSBarry Smith #else 419a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4205e42470aSBarry Smith #endif 4215e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4225e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4235e76c431SBarry Smith 424393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4255334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 42678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 427313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 428406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 429393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 43094a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 431ad922ac9SBarry Smith goto theend1; 4325e76c431SBarry Smith } 4335e76c431SBarry Smith 4345e76c431SBarry Smith /* Fit points with quadratic */ 435313b5538SBarry Smith lambda = 1.0; 436a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4375e76c431SBarry Smith lambdaprev = lambda; 4385e76c431SBarry Smith gnormprev = *gnorm; 439ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 440ddd12767SBarry Smith else lambda = lambdatemp; 441393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 442ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 443aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 444ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4455e42470aSBarry Smith #else 446ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4475e42470aSBarry Smith #endif 44878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 449cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 450406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 451393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 452ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 453ad922ac9SBarry Smith goto theend1; 4545e76c431SBarry Smith } 4555e76c431SBarry Smith 4565e76c431SBarry Smith /* Fit points with cubic */ 4575e76c431SBarry Smith count = 1; 4585e76c431SBarry Smith while (1) { 4595e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4602b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4612b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 462ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 463393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 464761aaf1bSLois Curfman McInnes *flag = -1; break; 4655e76c431SBarry Smith } 466406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 467406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 468ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4692b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4705e76c431SBarry Smith d = b*b - 3*a*initslope; 4715e76c431SBarry Smith if (d < 0.0) d = 0.0; 4725e76c431SBarry Smith if (a == 0.0) { 4735e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4745e76c431SBarry Smith } else { 4755e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4765e76c431SBarry Smith } 4775e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4785e76c431SBarry Smith lambdatemp = .5*lambda; 4795e76c431SBarry Smith } 4805e76c431SBarry Smith lambdaprev = lambda; 4815e76c431SBarry Smith gnormprev = *gnorm; 4825e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4835e76c431SBarry Smith lambda = .1*lambda; 4845e76c431SBarry Smith } 4855e76c431SBarry Smith else lambda = lambdatemp; 486393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w );CHKERRQ(ierr); 487ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 488aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 489ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 490393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4915e42470aSBarry Smith #else 492ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4935e42470aSBarry Smith #endif 49478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 495cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 496406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 497393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 498ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4992715565aSLois Curfman McInnes goto theend1; 5002b022350SLois Curfman McInnes } else { 5012b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5025e76c431SBarry Smith } 5035e76c431SBarry Smith count++; 5045e76c431SBarry Smith } 505ad922ac9SBarry Smith theend1: 5068f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5078f99978dSLois Curfman McInnes if (neP->CheckStep) { 5088f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5098f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5108f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5118f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5128f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5138f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5148f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5158f99978dSLois Curfman McInnes } 5168f99978dSLois Curfman McInnes } 5177857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5183a40ed3dSBarry Smith PetscFunctionReturn(0); 5195e76c431SBarry Smith } 52004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5215615d1e5SSatish Balay #undef __FUNC__ 5225615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 5234b828684SBarry Smith /*@C 524f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5255e76c431SBarry Smith 526c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 527c7afd0dbSLois Curfman McInnes 5285e42470aSBarry Smith Input Parameters: 529c7afd0dbSLois Curfman McInnes + snes - the SNES context 530af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5315e42470aSBarry Smith . x - current iterate 5325e42470aSBarry Smith . f - residual evaluated at x 5335e42470aSBarry Smith . y - search direction (contains new iterate on output) 5345e42470aSBarry Smith . w - work vector 535c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5365e42470aSBarry Smith 537c4a48953SLois Curfman McInnes Output Parameters: 538c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5395e42470aSBarry Smith . y - new iterate (contains search direction on input) 5405e42470aSBarry Smith . gnorm - 2-norm of g 5415e42470aSBarry Smith . ynorm - 2-norm of search length 542c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 543fee21e36SBarry Smith 544c4a48953SLois Curfman McInnes Options Database Key: 545c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 546c4a48953SLois Curfman McInnes 5475e42470aSBarry Smith Notes: 5486831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 5495e42470aSBarry Smith 55036851e7fSLois Curfman McInnes Level: advanced 55136851e7fSLois Curfman McInnes 552f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 553f59ffdedSLois Curfman McInnes 554af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5555e42470aSBarry Smith @*/ 556af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 557761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5585e76c431SBarry Smith { 559406556e6SLois Curfman McInnes /* 560406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 561406556e6SLois Curfman McInnes minimization problem: 562406556e6SLois Curfman McInnes min z(x): R^n -> R, 563406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 564406556e6SLois Curfman McInnes */ 56540011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 566aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5675e42470aSBarry Smith Scalar cinitslope,clambda; 5686b5873e3SBarry Smith #endif 5695e42470aSBarry Smith int ierr, count; 5706831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data; 5715334005bSBarry Smith Scalar mone = -1.0,scale; 5728f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5735e76c431SBarry Smith 5743a40ed3dSBarry Smith PetscFunctionBegin; 5757857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 576761aaf1bSLois Curfman McInnes *flag = 0; 5775e42470aSBarry Smith alpha = neP->alpha; 5785e42470aSBarry Smith maxstep = neP->maxstep; 5795e42470aSBarry Smith steptol = neP->steptol; 5805e76c431SBarry Smith 581cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 582b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 58394a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 584b37302c6SLois Curfman McInnes *gnorm = fnorm; 585b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 586b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 587ad922ac9SBarry Smith goto theend2; 58894a9d846SBarry Smith } 5895e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5905e42470aSBarry Smith scale = maxstep/(*ynorm); 591393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5925e42470aSBarry Smith *ynorm = maxstep; 5935e76c431SBarry Smith } 5945e42470aSBarry Smith minlambda = steptol/(*ynorm); 595a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 596aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 597a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 59892318cfeSSatish Balay initslope = PetscReal(cinitslope); 5995e42470aSBarry Smith #else 600a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6015e42470aSBarry Smith #endif 6025e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6035e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6045e42470aSBarry Smith 605393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6065334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 60778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 608cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 609406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 610393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 61194a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 612ad922ac9SBarry Smith goto theend2; 6135e42470aSBarry Smith } 6145e42470aSBarry Smith 6155e42470aSBarry Smith /* Fit points with quadratic */ 616313b5538SBarry Smith lambda = 1.0; 6175e42470aSBarry Smith count = 1; 6185e42470aSBarry Smith while (1) { 6195e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 620981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 621981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 622ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 623393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 624761aaf1bSLois Curfman McInnes *flag = -1; break; 6255e42470aSBarry Smith } 626a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 6275e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 6285e42470aSBarry Smith lambda = .1*lambda; 6295e42470aSBarry Smith } else lambda = lambdatemp; 630393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6315334005bSBarry Smith lambda = -lambda; 632aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 633393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6345e42470aSBarry Smith #else 635393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w);CHKERRQ(ierr); 6365e42470aSBarry Smith #endif 63778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 638cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 639406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 640393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 641981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 64206259719SBarry Smith break; 6435e42470aSBarry Smith } 6445e42470aSBarry Smith count++; 6455e42470aSBarry Smith } 646ad922ac9SBarry Smith theend2: 6478f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6488f99978dSLois Curfman McInnes if (neP->CheckStep) { 6498f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6508f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6518f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6528f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6538f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6548f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6558f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6568f99978dSLois Curfman McInnes } 6578f99978dSLois Curfman McInnes } 6587857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6593a40ed3dSBarry Smith PetscFunctionReturn(0); 6605e76c431SBarry Smith } 66104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6625615d1e5SSatish Balay #undef __FUNC__ 663d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 664c9e6a524SLois Curfman McInnes /*@C 66577c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 6666831982aSBarry Smith by the method SNESEQLS. 6675e76c431SBarry Smith 6685e76c431SBarry Smith Input Parameters: 669c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 670af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 671c7afd0dbSLois Curfman McInnes - func - pointer to int function 6725e76c431SBarry Smith 673fee21e36SBarry Smith Collective on SNES 674fee21e36SBarry Smith 675c4a48953SLois Curfman McInnes Available Routines: 676c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 677f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 678f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 679af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6805e76c431SBarry Smith 681c4a48953SLois Curfman McInnes Options Database Keys: 682af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 683c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 684c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 685c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 686c4a48953SLois Curfman McInnes 6875e76c431SBarry Smith Calling sequence of func: 688bd208895SLois Curfman McInnes .vb 689af2d14d2SLois Curfman McInnes func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 690bd208895SLois Curfman McInnes double fnorm, double *ynorm, double *gnorm, *flag) 691bd208895SLois Curfman McInnes .ve 6925e76c431SBarry Smith 6935e76c431SBarry Smith Input parameters for func: 694c7afd0dbSLois Curfman McInnes + snes - nonlinear context 695af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 6965e76c431SBarry Smith . x - current iterate 6975e76c431SBarry Smith . f - residual evaluated at x 6985e76c431SBarry Smith . y - search direction (contains new iterate on output) 6995e76c431SBarry Smith . w - work vector 700c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7015e76c431SBarry Smith 7025e76c431SBarry Smith Output parameters for func: 703c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7045e76c431SBarry Smith . y - new iterate (contains search direction on input) 7055e76c431SBarry Smith . gnorm - 2-norm of g 7065e76c431SBarry Smith . ynorm - 2-norm of search length 707c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 708761aaf1bSLois Curfman McInnes on failure. 709f59ffdedSLois Curfman McInnes 71036851e7fSLois Curfman McInnes Level: advanced 71136851e7fSLois Curfman McInnes 712f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 713f59ffdedSLois Curfman McInnes 7144ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7155e76c431SBarry Smith @*/ 716af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx) 7175e76c431SBarry Smith { 718af2d14d2SLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*); 71982bf6240SBarry Smith 7203a40ed3dSBarry Smith PetscFunctionBegin; 72194d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 72282bf6240SBarry Smith if (f) { 723af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 72482bf6240SBarry Smith } 7253a40ed3dSBarry Smith PetscFunctionReturn(0); 7265e76c431SBarry Smith } 72704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 728fb2e594dSBarry Smith EXTERN_C_BEGIN 72982bf6240SBarry Smith #undef __FUNC__ 73082bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 731af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 732af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 73382bf6240SBarry Smith { 73482bf6240SBarry Smith PetscFunctionBegin; 7356831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7366831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 73782bf6240SBarry Smith PetscFunctionReturn(0); 73882bf6240SBarry Smith } 739fb2e594dSBarry Smith EXTERN_C_END 74004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 741c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 742c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck" 743c8dd0c0dSLois Curfman McInnes /*@C 744530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7456831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 746c8dd0c0dSLois Curfman McInnes 747c8dd0c0dSLois Curfman McInnes Input Parameters: 748c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 749c8dd0c0dSLois Curfman McInnes . func - pointer to int function 750c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 751c8dd0c0dSLois Curfman McInnes 752c8dd0c0dSLois Curfman McInnes Collective on SNES 753c8dd0c0dSLois Curfman McInnes 754c8dd0c0dSLois Curfman McInnes Calling sequence of func: 755c8dd0c0dSLois Curfman McInnes .vb 756b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 757c8dd0c0dSLois Curfman McInnes .ve 758b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 759b1ae27eaSLois Curfman McInnes on failure. 760c8dd0c0dSLois Curfman McInnes 761c8dd0c0dSLois Curfman McInnes Input parameters for func: 762c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 763c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7648f99978dSLois Curfman McInnes - x - current candidate iterate 765c8dd0c0dSLois Curfman McInnes 766c8dd0c0dSLois Curfman McInnes Output parameters for func: 767c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 768c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 769c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 770c8dd0c0dSLois Curfman McInnes 771c8dd0c0dSLois Curfman McInnes Level: advanced 772c8dd0c0dSLois Curfman McInnes 7738f99978dSLois Curfman McInnes Notes: 774b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 775b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 776b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 777b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 778ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 7798f99978dSLois Curfman McInnes 780b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 781b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 782b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 783ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 784ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 785ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 786ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 787b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 7888f99978dSLois Curfman McInnes 789c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 790c8dd0c0dSLois Curfman McInnes 791c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 792c8dd0c0dSLois Curfman McInnes @*/ 7938f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 794c8dd0c0dSLois Curfman McInnes { 7958f99978dSLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 796c8dd0c0dSLois Curfman McInnes 797c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 798c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 799c8dd0c0dSLois Curfman McInnes if (f) { 800c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 801c8dd0c0dSLois Curfman McInnes } 802c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 803c8dd0c0dSLois Curfman McInnes } 804c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 805c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 806c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 807c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS" 8088f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 809c8dd0c0dSLois Curfman McInnes { 810c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8116831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8126831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 813c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 814c8dd0c0dSLois Curfman McInnes } 815c8dd0c0dSLois Curfman McInnes EXTERN_C_END 816c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 81704d965bbSLois Curfman McInnes /* 8186831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 81982bf6240SBarry Smith 82004d965bbSLois Curfman McInnes Input Parameter: 82104d965bbSLois Curfman McInnes . snes - the SNES context 82204d965bbSLois Curfman McInnes 82304d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 82404d965bbSLois Curfman McInnes */ 8255615d1e5SSatish Balay #undef __FUNC__ 826d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 827f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8285e42470aSBarry Smith { 8296831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 830d132466eSBarry Smith int ierr; 8316daaf66cSBarry Smith 8323a40ed3dSBarry Smith PetscFunctionBegin; 8336831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 834d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 835d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 836d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 837d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8383a40ed3dSBarry Smith PetscFunctionReturn(0); 8395e42470aSBarry Smith } 84004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 84104d965bbSLois Curfman McInnes /* 8426831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 84304d965bbSLois Curfman McInnes 84404d965bbSLois Curfman McInnes Input Parameters: 84504d965bbSLois Curfman McInnes . SNES - the SNES context 84604d965bbSLois Curfman McInnes . viewer - visualization context 84704d965bbSLois Curfman McInnes 84804d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 84904d965bbSLois Curfman McInnes */ 8505615d1e5SSatish Balay #undef __FUNC__ 851d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 852e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 853a935fc98SLois Curfman McInnes { 8546831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 85519bcc07fSBarry Smith char *cstr; 85651695f54SLois Curfman McInnes int ierr; 8576831982aSBarry Smith PetscTruth isascii; 858a935fc98SLois Curfman McInnes 8593a40ed3dSBarry Smith PetscFunctionBegin; 8606831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8610f5bd95cSBarry Smith if (isascii) { 86219bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 86319bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 86419bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 86519bcc07fSBarry Smith else cstr = "unknown"; 866d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 867d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 8685cd90555SBarry Smith } else { 8690f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 87019bcc07fSBarry Smith } 8713a40ed3dSBarry Smith PetscFunctionReturn(0); 872a935fc98SLois Curfman McInnes } 87304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 87404d965bbSLois Curfman McInnes /* 8756831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 87604d965bbSLois Curfman McInnes 87704d965bbSLois Curfman McInnes Input Parameter: 87804d965bbSLois Curfman McInnes . snes - the SNES context 87904d965bbSLois Curfman McInnes 88004d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 88104d965bbSLois Curfman McInnes */ 8825615d1e5SSatish Balay #undef __FUNC__ 8835615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 884f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8855e42470aSBarry Smith { 8866831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 8875e42470aSBarry Smith char ver[16]; 8885e42470aSBarry Smith double tmp; 889f1af5d2fSBarry Smith int ierr; 890f1af5d2fSBarry Smith PetscTruth flg; 8915e42470aSBarry Smith 8923a40ed3dSBarry Smith PetscFunctionBegin; 89309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr); 89425cdf11fSBarry Smith if (flg) { 8955e42470aSBarry Smith ls->alpha = tmp; 8965e42470aSBarry Smith } 89709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr); 89825cdf11fSBarry Smith if (flg) { 8995e42470aSBarry Smith ls->maxstep = tmp; 9005e42470aSBarry Smith } 90109d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr); 90225cdf11fSBarry Smith if (flg) { 9035e42470aSBarry Smith ls->steptol = tmp; 9045e42470aSBarry Smith } 90509d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg);CHKERRQ(ierr); 90625cdf11fSBarry Smith if (flg) { 907*c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9080f5bd95cSBarry Smith 909*c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 910*c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 911*c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 912*c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9130f5bd95cSBarry Smith 9140f5bd95cSBarry Smith if (isbasic) { 915af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9160f5bd95cSBarry Smith } else if (isnonorms) { 917af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9180f5bd95cSBarry Smith } else if (isquad) { 919af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9200f5bd95cSBarry Smith } else if (iscubic) { 921af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9225e42470aSBarry Smith } 923a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9245e42470aSBarry Smith } 9253a40ed3dSBarry Smith PetscFunctionReturn(0); 9265e42470aSBarry Smith } 92704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92804d965bbSLois Curfman McInnes /* 9296831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9306831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 93104d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 93204d965bbSLois Curfman McInnes 93304d965bbSLois Curfman McInnes 93404d965bbSLois Curfman McInnes Input Parameter: 93504d965bbSLois Curfman McInnes . snes - the SNES context 93604d965bbSLois Curfman McInnes 93704d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 93804d965bbSLois Curfman McInnes */ 939fb2e594dSBarry Smith EXTERN_C_BEGIN 9405615d1e5SSatish Balay #undef __FUNC__ 9415615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 942f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9435e42470aSBarry Smith { 94482bf6240SBarry Smith int ierr; 9456831982aSBarry Smith SNES_EQ_LS *neP; 9465e42470aSBarry Smith 9473a40ed3dSBarry Smith PetscFunctionBegin; 948a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 949a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 950a8c6a408SBarry Smith } 95182bf6240SBarry Smith 952f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 953f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 954f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 955f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 956f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 957f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 958f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9595baf8537SBarry Smith snes->nwork = 0; 9605e42470aSBarry Smith 9616831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 9626831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 9635e42470aSBarry Smith snes->data = (void *) neP; 9645e42470aSBarry Smith neP->alpha = 1.e-4; 9655e42470aSBarry Smith neP->maxstep = 1.e8; 9665e42470aSBarry Smith neP->steptol = 1.e-12; 9675e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 968c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 969c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 970c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 97182bf6240SBarry Smith 972f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 973e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 974f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 975c8dd0c0dSLois Curfman McInnes (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 97682bf6240SBarry Smith 9773a40ed3dSBarry Smith PetscFunctionReturn(0); 9785e42470aSBarry Smith } 979fb2e594dSBarry Smith EXTERN_C_END 98082bf6240SBarry Smith 98182bf6240SBarry Smith 98282bf6240SBarry Smith 983