1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*8f99978dSLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.126 1999/03/14 17:54:40 curfman Exp curfman $"; 35e76c431SBarry Smith #endif 45e76c431SBarry Smith 570f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 65e42470aSBarry Smith 704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 804d965bbSLois Curfman McInnes 904d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 1004d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 1104d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 1204d965bbSLois Curfman McInnes respectively. 1304d965bbSLois Curfman McInnes 14fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 1504d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 1604d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 17fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 1804d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 1904d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 2004d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 2104d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 2204d965bbSLois Curfman McInnes These routines are actually called via the common user interface 2304d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 2404d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 2504d965bbSLois Curfman McInnes for all nonlinear solvers. 2604d965bbSLois Curfman McInnes 2704d965bbSLois Curfman McInnes Another key routine is: 2804d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 2904d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 3004d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 3104d965bbSLois Curfman McInnes SNESSolve() if necessary. 3204d965bbSLois Curfman McInnes 3304d965bbSLois Curfman McInnes Additional basic routines are: 3404d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 3504d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 3604d965bbSLois Curfman McInnes have actually been used. 3704d965bbSLois Curfman McInnes These are called by application codes via the interface routines 3804d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 3904d965bbSLois Curfman McInnes 4004d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 4104d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 4204d965bbSLois Curfman McInnes above description applies to these categories also. 4304d965bbSLois Curfman McInnes 4404d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 455e42470aSBarry Smith /* 4604d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 4704d965bbSLois Curfman McInnes method with a line search. 485e76c431SBarry Smith 4904d965bbSLois Curfman McInnes Input Parameters: 5004d965bbSLois Curfman McInnes . snes - the SNES context 515e76c431SBarry Smith 5204d965bbSLois Curfman McInnes Output Parameter: 53c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 5404d965bbSLois Curfman McInnes 5504d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 565e76c431SBarry Smith 575e76c431SBarry Smith Notes: 585e76c431SBarry Smith This implements essentially a truncated Newton method with a 595e76c431SBarry Smith line search. By default a cubic backtracking line search 605e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 615e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 62393d2d9aSLois Curfman McInnes and Schnabel. 635e76c431SBarry Smith */ 645615d1e5SSatish Balay #undef __FUNC__ 655615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS" 66f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 675e76c431SBarry Smith { 685e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 6984c569c9SBarry Smith int maxits, i, ierr, lits, lsfail; 70112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 7184c569c9SBarry Smith double fnorm, gnorm, xnorm, ynorm; 725e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 735e76c431SBarry Smith 743a40ed3dSBarry Smith PetscFunctionBegin; 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|| */ 8406d0569aSBarry Smith PetscAMSTakeAccess(snes); 8584c569c9SBarry Smith snes->iter = 0; 865e42470aSBarry Smith snes->norm = fnorm; 8706d0569aSBarry Smith PetscAMSGrantAccess(snes); 8884c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 8994a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 905e76c431SBarry Smith 913a40ed3dSBarry Smith if (fnorm < snes->atol) {*outits = 0; 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); 1027a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(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); 112761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 1135e76c431SBarry Smith 11439e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 11539e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1165e76c431SBarry Smith fnorm = gnorm; 1175e76c431SBarry Smith 11806d0569aSBarry Smith PetscAMSTakeAccess(snes); 11984c569c9SBarry Smith snes->iter = i+1; 1205e42470aSBarry Smith snes->norm = fnorm; 12106d0569aSBarry Smith PetscAMSGrantAccess(snes); 12284c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 12394a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1245e76c431SBarry Smith 1255e76c431SBarry Smith /* Test for convergence */ 12629e0b56fSBarry Smith if (snes->converged) { 12729e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr); /* xnorm = || X || */ 128bbb6d6a8SBarry Smith if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) { 12916c95cacSBarry Smith break; 13016c95cacSBarry Smith } 13116c95cacSBarry Smith } 13229e0b56fSBarry Smith } 13339e2f89bSBarry Smith if (X != snes->vec_sol) { 134393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr); 13539e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 13639e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 13739e2f89bSBarry Smith } 13852392280SLois Curfman McInnes if (i == maxits) { 139981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 14052392280SLois Curfman McInnes i--; 14152392280SLois Curfman McInnes } 1425e42470aSBarry Smith *outits = i+1; 1433a40ed3dSBarry Smith PetscFunctionReturn(0); 1445e76c431SBarry Smith } 14504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 14604d965bbSLois Curfman McInnes /* 14704d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 14804d965bbSLois Curfman McInnes of the SNES_EQ_LS nonlinear solver. 14904d965bbSLois Curfman McInnes 15004d965bbSLois Curfman McInnes Input Parameter: 15104d965bbSLois Curfman McInnes . snes - the SNES context 15204d965bbSLois Curfman McInnes . x - the solution vector 15304d965bbSLois Curfman McInnes 15404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 15504d965bbSLois Curfman McInnes 15604d965bbSLois Curfman McInnes Notes: 15704d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 15804d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 15904d965bbSLois Curfman McInnes the call to SNESSolve(). 16004d965bbSLois Curfman McInnes */ 1615615d1e5SSatish Balay #undef __FUNC__ 1625615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS" 163f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1645e76c431SBarry Smith { 1655e42470aSBarry Smith int ierr; 1663a40ed3dSBarry Smith 1673a40ed3dSBarry Smith PetscFunctionBegin; 16881b6cf68SLois Curfman McInnes snes->nwork = 4; 169d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1705e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 17181b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1723a40ed3dSBarry Smith PetscFunctionReturn(0); 1735e76c431SBarry Smith } 17404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 17504d965bbSLois Curfman McInnes /* 17604d965bbSLois Curfman McInnes SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created 17704d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 17804d965bbSLois Curfman McInnes 17904d965bbSLois Curfman McInnes Input Parameter: 18004d965bbSLois Curfman McInnes . snes - the SNES context 18104d965bbSLois Curfman McInnes 182de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 18304d965bbSLois Curfman McInnes */ 1845615d1e5SSatish Balay #undef __FUNC__ 185d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS" 186e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1875e76c431SBarry Smith { 188393d2d9aSLois Curfman McInnes int ierr; 1893a40ed3dSBarry Smith 1903a40ed3dSBarry Smith PetscFunctionBegin; 1915baf8537SBarry Smith if (snes->nwork) { 1924b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr); 19321c89e3eSBarry Smith } 1940452661fSBarry Smith PetscFree(snes->data); 1953a40ed3dSBarry Smith PetscFunctionReturn(0); 1965e76c431SBarry Smith } 19704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1985615d1e5SSatish Balay #undef __FUNC__ 1995615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch" 20082bf6240SBarry Smith 2014b828684SBarry Smith /*@C 2025e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2035e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2045e76c431SBarry Smith to serve as a template and is not recommended for general use. 2055e76c431SBarry Smith 206c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 207c7afd0dbSLois Curfman McInnes 2085e76c431SBarry Smith Input Parameters: 209c7afd0dbSLois Curfman McInnes + snes - nonlinear context 210af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2115e76c431SBarry Smith . x - current iterate 2125e76c431SBarry Smith . f - residual evaluated at x 2135e76c431SBarry Smith . y - search direction (contains new iterate on output) 2145e76c431SBarry Smith . w - work vector 215c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2165e76c431SBarry Smith 217c4a48953SLois Curfman McInnes Output Parameters: 218c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2195e76c431SBarry Smith . y - new iterate (contains search direction on input) 2205e76c431SBarry Smith . gnorm - 2-norm of g 2215e76c431SBarry Smith . ynorm - 2-norm of search length 222c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 223fee21e36SBarry Smith 224c4a48953SLois Curfman McInnes Options Database Key: 225c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 226c4a48953SLois Curfman McInnes 22736851e7fSLois Curfman McInnes Level: advanced 22836851e7fSLois Curfman McInnes 22928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 23028ae5a21SLois Curfman McInnes 231f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 232af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2335e76c431SBarry Smith @*/ 234af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 235761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 2365e76c431SBarry Smith { 2375e42470aSBarry Smith int ierr; 2385334005bSBarry Smith Scalar mone = -1.0; 239*8f99978dSLois Curfman McInnes SNES_LS *neP = (SNES_LS *) snes->data; 240*8f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2415334005bSBarry Smith 2423a40ed3dSBarry Smith PetscFunctionBegin; 243761aaf1bSLois Curfman McInnes *flag = 0; 2447857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 245ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 246*8f99978dSLois Curfman McInnes if (neP->CheckStep) { 247*8f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 248*8f99978dSLois Curfman McInnes } 249ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 250*8f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 251*8f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 252*8f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 253*8f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr); 2547857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 2565e76c431SBarry Smith } 25704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2585615d1e5SSatish Balay #undef __FUNC__ 25929e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 26082bf6240SBarry Smith 26129e0b56fSBarry Smith /*@C 26229e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 26329e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 26429e0b56fSBarry Smith even compute the norm of the function or search direction; this 26529e0b56fSBarry Smith is intended only when you know the full step is fine and are 26629e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 267c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 268c7afd0dbSLois Curfman McInnes 269c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 27029e0b56fSBarry Smith 27129e0b56fSBarry Smith Input Parameters: 272c7afd0dbSLois Curfman McInnes + snes - nonlinear context 273af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 27429e0b56fSBarry Smith . x - current iterate 27529e0b56fSBarry Smith . f - residual evaluated at x 27629e0b56fSBarry Smith . y - search direction (contains new iterate on output) 27729e0b56fSBarry Smith . w - work vector 278c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 27929e0b56fSBarry Smith 28029e0b56fSBarry Smith Output Parameters: 281c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 28229e0b56fSBarry Smith . gnorm - not changed 28329e0b56fSBarry Smith . ynorm - not changed 284c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 285fee21e36SBarry Smith 28629e0b56fSBarry Smith Options Database Key: 287c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 28829e0b56fSBarry Smith 28936851e7fSLois Curfman McInnes Level: advanced 29036851e7fSLois Curfman McInnes 29129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 29229e0b56fSBarry Smith 29329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 29429e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 29529e0b56fSBarry Smith @*/ 296af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 29729e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 29829e0b56fSBarry Smith { 29929e0b56fSBarry Smith int ierr; 30029e0b56fSBarry Smith Scalar mone = -1.0; 301*8f99978dSLois Curfman McInnes SNES_LS *neP = (SNES_LS *) snes->data; 302*8f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 30329e0b56fSBarry Smith 3043a40ed3dSBarry Smith PetscFunctionBegin; 30529e0b56fSBarry Smith *flag = 0; 30629e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 30729e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 308*8f99978dSLois Curfman McInnes if (neP->CheckStep) { 309*8f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 310*8f99978dSLois Curfman McInnes } 31129e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 31229e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 31429e0b56fSBarry Smith } 31504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 31629e0b56fSBarry Smith #undef __FUNC__ 3175615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3184b828684SBarry Smith /*@C 319f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3205e76c431SBarry Smith 321c7afd0dbSLois Curfman McInnes Collective on SNES 322c7afd0dbSLois Curfman McInnes 3235e76c431SBarry Smith Input Parameters: 324c7afd0dbSLois Curfman McInnes + snes - nonlinear context 325af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3265e76c431SBarry Smith . x - current iterate 3275e76c431SBarry Smith . f - residual evaluated at x 3285e76c431SBarry Smith . y - search direction (contains new iterate on output) 3295e76c431SBarry Smith . w - work vector 330c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3315e76c431SBarry Smith 332393d2d9aSLois Curfman McInnes Output Parameters: 333c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3345e76c431SBarry Smith . y - new iterate (contains search direction on input) 3355e76c431SBarry Smith . gnorm - 2-norm of g 3365e76c431SBarry Smith . ynorm - 2-norm of search length 337c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 338fee21e36SBarry Smith 339c4a48953SLois Curfman McInnes Options Database Key: 340c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 341c4a48953SLois Curfman McInnes 3425e76c431SBarry Smith Notes: 3435e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3445e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3455e76c431SBarry Smith 34636851e7fSLois Curfman McInnes Level: advanced 34736851e7fSLois Curfman McInnes 34828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 34928ae5a21SLois Curfman McInnes 350af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3515e76c431SBarry Smith @*/ 352af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 353761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3545e76c431SBarry Smith { 355406556e6SLois Curfman McInnes /* 356406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 357406556e6SLois Curfman McInnes minimization problem: 358406556e6SLois Curfman McInnes min z(x): R^n -> R, 359406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 360406556e6SLois Curfman McInnes */ 361406556e6SLois Curfman McInnes 362ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 363ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 3643a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 3655e42470aSBarry Smith Scalar cinitslope, clambda; 3666b5873e3SBarry Smith #endif 3675e42470aSBarry Smith int ierr, count; 3685e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3695334005bSBarry Smith Scalar mone = -1.0, scale; 370*8f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3715e76c431SBarry Smith 3723a40ed3dSBarry Smith PetscFunctionBegin; 3737857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 374761aaf1bSLois Curfman McInnes *flag = 0; 3755e76c431SBarry Smith alpha = neP->alpha; 3765e76c431SBarry Smith maxstep = neP->maxstep; 3775e76c431SBarry Smith steptol = neP->steptol; 3785e76c431SBarry Smith 379cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); 380a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 381a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 382a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 383a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 384a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 385ad922ac9SBarry Smith goto theend1; 38694a9d846SBarry Smith } 3875e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3885e42470aSBarry Smith scale = maxstep/(*ynorm); 3893a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 39092318cfeSSatish Balay PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale)); 3916b5873e3SBarry Smith #else 39294a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 3936b5873e3SBarry Smith #endif 394393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 3955e76c431SBarry Smith *ynorm = maxstep; 3965e76c431SBarry Smith } 3975e76c431SBarry Smith minlambda = steptol/(*ynorm); 398a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 3993a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 400a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 40192318cfeSSatish Balay initslope = PetscReal(cinitslope); 4025e42470aSBarry Smith #else 403a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 4045e42470aSBarry Smith #endif 4055e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4065e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4075e76c431SBarry Smith 408393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 4095334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 41078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 411cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); 412406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 413393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 41494a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 415ad922ac9SBarry Smith goto theend1; 4165e76c431SBarry Smith } 4175e76c431SBarry Smith 4185e76c431SBarry Smith /* Fit points with quadratic */ 4195e76c431SBarry Smith lambda = 1.0; count = 0; 420a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4215e76c431SBarry Smith lambdaprev = lambda; 4225e76c431SBarry Smith gnormprev = *gnorm; 423ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 424ddd12767SBarry Smith else lambda = lambdatemp; 425393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 426ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4273a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 428ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4295e42470aSBarry Smith #else 430ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4315e42470aSBarry Smith #endif 43278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 433cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 434406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 435393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 436ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 437ad922ac9SBarry Smith goto theend1; 4385e76c431SBarry Smith } 4395e76c431SBarry Smith 4405e76c431SBarry Smith /* Fit points with cubic */ 4415e76c431SBarry Smith count = 1; 4425e76c431SBarry Smith while (1) { 4435e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4442b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4452b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 446ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 447393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 448761aaf1bSLois Curfman McInnes *flag = -1; break; 4495e76c431SBarry Smith } 450406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 451406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 452ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4532b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4545e76c431SBarry Smith d = b*b - 3*a*initslope; 4555e76c431SBarry Smith if (d < 0.0) d = 0.0; 4565e76c431SBarry Smith if (a == 0.0) { 4575e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4585e76c431SBarry Smith } else { 4595e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4605e76c431SBarry Smith } 4615e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4625e76c431SBarry Smith lambdatemp = .5*lambda; 4635e76c431SBarry Smith } 4645e76c431SBarry Smith lambdaprev = lambda; 4655e76c431SBarry Smith gnormprev = *gnorm; 4665e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4675e76c431SBarry Smith lambda = .1*lambda; 4685e76c431SBarry Smith } 4695e76c431SBarry Smith else lambda = lambdatemp; 470393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 471ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4723a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 473ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 474393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4755e42470aSBarry Smith #else 476ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4775e42470aSBarry Smith #endif 47878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 479cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 480406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 481393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 482ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4832715565aSLois Curfman McInnes goto theend1; 4842b022350SLois Curfman McInnes } else { 4852b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 4865e76c431SBarry Smith } 4875e76c431SBarry Smith count++; 4885e76c431SBarry Smith } 489ad922ac9SBarry Smith theend1: 490*8f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 491*8f99978dSLois Curfman McInnes if (neP->CheckStep) { 492*8f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 493*8f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 494*8f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); 495*8f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); 496*8f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); 497*8f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 498*8f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr); 499*8f99978dSLois Curfman McInnes } 500*8f99978dSLois Curfman McInnes } 5017857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5023a40ed3dSBarry Smith PetscFunctionReturn(0); 5035e76c431SBarry Smith } 50404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5055615d1e5SSatish Balay #undef __FUNC__ 5065615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 5074b828684SBarry Smith /*@C 508f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5095e76c431SBarry Smith 510c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 511c7afd0dbSLois Curfman McInnes 5125e42470aSBarry Smith Input Parameters: 513c7afd0dbSLois Curfman McInnes + snes - the SNES context 514af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5155e42470aSBarry Smith . x - current iterate 5165e42470aSBarry Smith . f - residual evaluated at x 5175e42470aSBarry Smith . y - search direction (contains new iterate on output) 5185e42470aSBarry Smith . w - work vector 519c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5205e42470aSBarry Smith 521c4a48953SLois Curfman McInnes Output Parameters: 522c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5235e42470aSBarry Smith . y - new iterate (contains search direction on input) 5245e42470aSBarry Smith . gnorm - 2-norm of g 5255e42470aSBarry Smith . ynorm - 2-norm of search length 526c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 527fee21e36SBarry Smith 528c4a48953SLois Curfman McInnes Options Database Key: 529c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 530c4a48953SLois Curfman McInnes 5315e42470aSBarry Smith Notes: 532c7afd0dbSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method. 5335e42470aSBarry Smith 53436851e7fSLois Curfman McInnes Level: advanced 53536851e7fSLois Curfman McInnes 536f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 537f59ffdedSLois Curfman McInnes 538af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5395e42470aSBarry Smith @*/ 540af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 541761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5425e76c431SBarry Smith { 543406556e6SLois Curfman McInnes /* 544406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 545406556e6SLois Curfman McInnes minimization problem: 546406556e6SLois Curfman McInnes min z(x): R^n -> R, 547406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 548406556e6SLois Curfman McInnes */ 54940011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 5503a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 5515e42470aSBarry Smith Scalar cinitslope,clambda; 5526b5873e3SBarry Smith #endif 5535e42470aSBarry Smith int ierr, count; 5545e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 5555334005bSBarry Smith Scalar mone = -1.0,scale; 556*8f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5575e76c431SBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 5597857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 560761aaf1bSLois Curfman McInnes *flag = 0; 5615e42470aSBarry Smith alpha = neP->alpha; 5625e42470aSBarry Smith maxstep = neP->maxstep; 5635e42470aSBarry Smith steptol = neP->steptol; 5645e76c431SBarry Smith 565cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 566b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 56794a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 568b37302c6SLois Curfman McInnes *gnorm = fnorm; 569b37302c6SLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 570b37302c6SLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 571ad922ac9SBarry Smith goto theend2; 57294a9d846SBarry Smith } 5735e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5745e42470aSBarry Smith scale = maxstep/(*ynorm); 575393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 5765e42470aSBarry Smith *ynorm = maxstep; 5775e76c431SBarry Smith } 5785e42470aSBarry Smith minlambda = steptol/(*ynorm); 579a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 5803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 581a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 58292318cfeSSatish Balay initslope = PetscReal(cinitslope); 5835e42470aSBarry Smith #else 584a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 5855e42470aSBarry Smith #endif 5865e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 5875e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 5885e42470aSBarry Smith 589393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 5905334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 59178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 592cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 593406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 594393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 59594a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 596ad922ac9SBarry Smith goto theend2; 5975e42470aSBarry Smith } 5985e42470aSBarry Smith 5995e42470aSBarry Smith /* Fit points with quadratic */ 6005e42470aSBarry Smith lambda = 1.0; count = 0; 6015e42470aSBarry Smith count = 1; 6025e42470aSBarry Smith while (1) { 6035e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 604981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 605981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 606ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 607393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 608761aaf1bSLois Curfman McInnes *flag = -1; break; 6095e42470aSBarry Smith } 610a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 6115e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 6125e42470aSBarry Smith lambda = .1*lambda; 6135e42470aSBarry Smith } else lambda = lambdatemp; 614393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 6155334005bSBarry Smith lambda = -lambda; 6163a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 617393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 6185e42470aSBarry Smith #else 619393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 6205e42470aSBarry Smith #endif 62178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 622cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 623406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 624393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 625981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 62606259719SBarry Smith break; 6275e42470aSBarry Smith } 6285e42470aSBarry Smith count++; 6295e42470aSBarry Smith } 630ad922ac9SBarry Smith theend2: 631*8f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 632*8f99978dSLois Curfman McInnes if (neP->CheckStep) { 633*8f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 634*8f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 635*8f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); 636*8f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); 637*8f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); 638*8f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 639*8f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr); 640*8f99978dSLois Curfman McInnes } 641*8f99978dSLois Curfman McInnes } 6427857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6433a40ed3dSBarry Smith PetscFunctionReturn(0); 6445e76c431SBarry Smith } 64504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6465615d1e5SSatish Balay #undef __FUNC__ 647d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 648c9e6a524SLois Curfman McInnes /*@C 64977c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 650f63b844aSLois Curfman McInnes by the method SNES_EQ_LS. 6515e76c431SBarry Smith 6525e76c431SBarry Smith Input Parameters: 653c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 654af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 655c7afd0dbSLois Curfman McInnes - func - pointer to int function 6565e76c431SBarry Smith 657fee21e36SBarry Smith Collective on SNES 658fee21e36SBarry Smith 659c4a48953SLois Curfman McInnes Available Routines: 660c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 661f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 662f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 663af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6645e76c431SBarry Smith 665c4a48953SLois Curfman McInnes Options Database Keys: 666af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 667c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 668c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 669c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 670c4a48953SLois Curfman McInnes 6715e76c431SBarry Smith Calling sequence of func: 672bd208895SLois Curfman McInnes .vb 673af2d14d2SLois Curfman McInnes func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 674bd208895SLois Curfman McInnes double fnorm, double *ynorm, double *gnorm, *flag) 675bd208895SLois Curfman McInnes .ve 6765e76c431SBarry Smith 6775e76c431SBarry Smith Input parameters for func: 678c7afd0dbSLois Curfman McInnes + snes - nonlinear context 679af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 6805e76c431SBarry Smith . x - current iterate 6815e76c431SBarry Smith . f - residual evaluated at x 6825e76c431SBarry Smith . y - search direction (contains new iterate on output) 6835e76c431SBarry Smith . w - work vector 684c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6855e76c431SBarry Smith 6865e76c431SBarry Smith Output parameters for func: 687c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6885e76c431SBarry Smith . y - new iterate (contains search direction on input) 6895e76c431SBarry Smith . gnorm - 2-norm of g 6905e76c431SBarry Smith . ynorm - 2-norm of search length 691c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 692761aaf1bSLois Curfman McInnes on failure. 693f59ffdedSLois Curfman McInnes 69436851e7fSLois Curfman McInnes Level: advanced 69536851e7fSLois Curfman McInnes 696f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 697f59ffdedSLois Curfman McInnes 698*8f99978dSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck() 6995e76c431SBarry Smith @*/ 700af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx) 7015e76c431SBarry Smith { 702af2d14d2SLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*); 70382bf6240SBarry Smith 7043a40ed3dSBarry Smith PetscFunctionBegin; 70594d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 70682bf6240SBarry Smith if (f) { 707af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 70882bf6240SBarry Smith } 7093a40ed3dSBarry Smith PetscFunctionReturn(0); 7105e76c431SBarry Smith } 71104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 712fb2e594dSBarry Smith EXTERN_C_BEGIN 71382bf6240SBarry Smith #undef __FUNC__ 71482bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 715af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 716af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 71782bf6240SBarry Smith { 71882bf6240SBarry Smith PetscFunctionBegin; 71982bf6240SBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 720af2d14d2SLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 72182bf6240SBarry Smith PetscFunctionReturn(0); 72282bf6240SBarry Smith } 723fb2e594dSBarry Smith EXTERN_C_END 72404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 725c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 726c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck" 727c8dd0c0dSLois Curfman McInnes /*@C 728*8f99978dSLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the new iterate computed 729*8f99978dSLois Curfman McInnes by the line search routine in the Newton-based method SNES_EQ_LS. 730c8dd0c0dSLois Curfman McInnes 731c8dd0c0dSLois Curfman McInnes Input Parameters: 732c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 733c8dd0c0dSLois Curfman McInnes . func - pointer to int function 734c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 735c8dd0c0dSLois Curfman McInnes 736c8dd0c0dSLois Curfman McInnes Collective on SNES 737c8dd0c0dSLois Curfman McInnes 738c8dd0c0dSLois Curfman McInnes Calling sequence of func: 739c8dd0c0dSLois Curfman McInnes .vb 740*8f99978dSLois Curfman McInnes func (SNES snes, void *lsctx, Vec x, PetscTruth *flag) 741c8dd0c0dSLois Curfman McInnes .ve 742c8dd0c0dSLois Curfman McInnes 743c8dd0c0dSLois Curfman McInnes Input parameters for func: 744c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 745c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 746*8f99978dSLois Curfman McInnes - x - current candidate iterate 747c8dd0c0dSLois Curfman McInnes 748c8dd0c0dSLois Curfman McInnes Output parameters for func: 749c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 750c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 751c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 752c8dd0c0dSLois Curfman McInnes 753c8dd0c0dSLois Curfman McInnes Level: advanced 754c8dd0c0dSLois Curfman McInnes 755*8f99978dSLois Curfman McInnes Notes: 756*8f99978dSLois Curfman McInnes The user-defined line search checking routine is available for 757*8f99978dSLois Curfman McInnes use in conjunction with SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 758*8f99978dSLois Curfman McInnes which (1) compute a candidate iterate u_{i+1}, (2) pass control to the 759*8f99978dSLois Curfman McInnes checking routine, and then (3) compute the corresponding function f(u_{i+1}) 760*8f99978dSLois Curfman McInnes with the (possibly altered) iterate u_{i+1}. 761*8f99978dSLois Curfman McInnes 762*8f99978dSLois Curfman McInnes The user-defined line search checking routine is also available for 763*8f99978dSLois Curfman McInnes use in conjunction with SNESQuadraticLineSearch() and SNESCubicLineSearch(). 764*8f99978dSLois Curfman McInnes These routines (1) compute a candidate iterate u_{i+1} as well as a 765*8f99978dSLois Curfman McInnes candidate function f(u_{i+1}), (2) pass control to the checking routine, 766*8f99978dSLois Curfman McInnes and then (3) force a re-evaluation of f(u_{i+1}) if any changes were made 767*8f99978dSLois Curfman McInnes to the candidate iterate in the checking routine (as indicated by 768*8f99978dSLois Curfman McInnes flag=PETSC_TRUE). The overhead of this re-evaluation can be costly, so 769*8f99978dSLois Curfman McInnes this feature with caution. 770*8f99978dSLois Curfman McInnes 771c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 772c8dd0c0dSLois Curfman McInnes 773c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 774c8dd0c0dSLois Curfman McInnes @*/ 775*8f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 776c8dd0c0dSLois Curfman McInnes { 777*8f99978dSLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 778c8dd0c0dSLois Curfman McInnes 779c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 780c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 781c8dd0c0dSLois Curfman McInnes if (f) { 782c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 783c8dd0c0dSLois Curfman McInnes } 784c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 785c8dd0c0dSLois Curfman McInnes } 786c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 787c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 788c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 789c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS" 790*8f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 791c8dd0c0dSLois Curfman McInnes { 792c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 793c8dd0c0dSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 794c8dd0c0dSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 795c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 796c8dd0c0dSLois Curfman McInnes } 797c8dd0c0dSLois Curfman McInnes EXTERN_C_END 798c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 79904d965bbSLois Curfman McInnes /* 80004d965bbSLois Curfman McInnes SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 80182bf6240SBarry Smith 80204d965bbSLois Curfman McInnes Input Parameter: 80304d965bbSLois Curfman McInnes . snes - the SNES context 80404d965bbSLois Curfman McInnes 80504d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 80604d965bbSLois Curfman McInnes */ 8075615d1e5SSatish Balay #undef __FUNC__ 808d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 809f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8105e42470aSBarry Smith { 8115e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 8126daaf66cSBarry Smith 8133a40ed3dSBarry Smith PetscFunctionBegin; 81476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n"); 815af2d14d2SLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 81676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 81776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 81876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol); 8193a40ed3dSBarry Smith PetscFunctionReturn(0); 8205e42470aSBarry Smith } 82104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 82204d965bbSLois Curfman McInnes /* 823de49b36dSLois Curfman McInnes SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure. 82404d965bbSLois Curfman McInnes 82504d965bbSLois Curfman McInnes Input Parameters: 82604d965bbSLois Curfman McInnes . SNES - the SNES context 82704d965bbSLois Curfman McInnes . viewer - visualization context 82804d965bbSLois Curfman McInnes 82904d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 83004d965bbSLois Curfman McInnes */ 8315615d1e5SSatish Balay #undef __FUNC__ 832d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 833e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 834a935fc98SLois Curfman McInnes { 835a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 83619bcc07fSBarry Smith char *cstr; 83751695f54SLois Curfman McInnes int ierr; 83819bcc07fSBarry Smith ViewerType vtype; 839a935fc98SLois Curfman McInnes 8403a40ed3dSBarry Smith PetscFunctionBegin; 84119bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 8423f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 84319bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 84419bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 84519bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 84619bcc07fSBarry Smith else cstr = "unknown"; 8470ef38995SBarry Smith ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr); 8480ef38995SBarry Smith ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol); 8495cd90555SBarry Smith } else { 8505cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 85119bcc07fSBarry Smith } 8523a40ed3dSBarry Smith PetscFunctionReturn(0); 853a935fc98SLois Curfman McInnes } 85404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 85504d965bbSLois Curfman McInnes /* 85604d965bbSLois Curfman McInnes SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method. 85704d965bbSLois Curfman McInnes 85804d965bbSLois Curfman McInnes Input Parameter: 85904d965bbSLois Curfman McInnes . snes - the SNES context 86004d965bbSLois Curfman McInnes 86104d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 86204d965bbSLois Curfman McInnes */ 8635615d1e5SSatish Balay #undef __FUNC__ 8645615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 865f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8665e42470aSBarry Smith { 8675e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 8685e42470aSBarry Smith char ver[16]; 8695e42470aSBarry Smith double tmp; 87025cdf11fSBarry Smith int ierr,flg; 8715e42470aSBarry Smith 8723a40ed3dSBarry Smith PetscFunctionBegin; 87309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr); 87425cdf11fSBarry Smith if (flg) { 8755e42470aSBarry Smith ls->alpha = tmp; 8765e42470aSBarry Smith } 87709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr); 87825cdf11fSBarry Smith if (flg) { 8795e42470aSBarry Smith ls->maxstep = tmp; 8805e42470aSBarry Smith } 88109d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr); 88225cdf11fSBarry Smith if (flg) { 8835e42470aSBarry Smith ls->steptol = tmp; 8845e42470aSBarry Smith } 88509d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr); 88625cdf11fSBarry Smith if (flg) { 88748d91487SBarry Smith if (!PetscStrcmp(ver,"basic")) { 888af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 889b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"basicnonorms")) { 890af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 891b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"quadratic")) { 892af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 893b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"cubic")) { 894af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 8955e42470aSBarry Smith } 896a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 8975e42470aSBarry Smith } 8983a40ed3dSBarry Smith PetscFunctionReturn(0); 8995e42470aSBarry Smith } 90004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 90104d965bbSLois Curfman McInnes /* 90204d965bbSLois Curfman McInnes SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method, 90304d965bbSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 90404d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 90504d965bbSLois Curfman McInnes 90604d965bbSLois Curfman McInnes 90704d965bbSLois Curfman McInnes Input Parameter: 90804d965bbSLois Curfman McInnes . snes - the SNES context 90904d965bbSLois Curfman McInnes 91004d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 91104d965bbSLois Curfman McInnes */ 912fb2e594dSBarry Smith EXTERN_C_BEGIN 9135615d1e5SSatish Balay #undef __FUNC__ 9145615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 915f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9165e42470aSBarry Smith { 91782bf6240SBarry Smith int ierr; 9185e42470aSBarry Smith SNES_LS *neP; 9195e42470aSBarry Smith 9203a40ed3dSBarry Smith PetscFunctionBegin; 921a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 922a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 923a8c6a408SBarry Smith } 92482bf6240SBarry Smith 925f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 926f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 927f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 928f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 929f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 930f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 931f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9325baf8537SBarry Smith snes->nwork = 0; 9335e42470aSBarry Smith 9340452661fSBarry Smith neP = PetscNew(SNES_LS); CHKPTRQ(neP); 935ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 9365e42470aSBarry Smith snes->data = (void *) neP; 9375e42470aSBarry Smith neP->alpha = 1.e-4; 9385e42470aSBarry Smith neP->maxstep = 1.e8; 9395e42470aSBarry Smith neP->steptol = 1.e-12; 9405e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 941c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 942c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 943c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 94482bf6240SBarry Smith 94594d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 946e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 947c8dd0c0dSLois Curfman McInnes ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 948c8dd0c0dSLois Curfman McInnes (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 94982bf6240SBarry Smith 9503a40ed3dSBarry Smith PetscFunctionReturn(0); 9515e42470aSBarry Smith } 952fb2e594dSBarry Smith EXTERN_C_END 95382bf6240SBarry Smith 95482bf6240SBarry Smith 95582bf6240SBarry Smith 956