1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*8cbba510SBarry Smith static char vcid[] = "$Id: ls.c,v 1.128 1999/03/15 01:31:50 bsmith Exp bsmith $"; 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; 2398f99978dSLois Curfman McInnes SNES_LS *neP = (SNES_LS *) snes->data; 2408f99978dSLois 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 */ 2468f99978dSLois Curfman McInnes if (neP->CheckStep) { 2478f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 2488f99978dSLois Curfman McInnes } 249ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 2508f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 2518f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 2528f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 2538f99978dSLois 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 291*8cbba510SBarry Smith Notes: 292*8cbba510SBarry Smith You must run this with -snes_no_convergence_test or your own 293*8cbba510SBarry Smith custom test and use -snes_max_it nk or the SNES solver will generate 294*8cbba510SBarry Smith an error. 295*8cbba510SBarry Smith 296*8cbba510SBarry Smith Residual norm output printed from -snes_monitor will not be correct, 297*8cbba510SBarry Smith since they are not computed. 298*8cbba510SBarry Smith 29929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 30029e0b56fSBarry Smith 30129e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 30229e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 30329e0b56fSBarry Smith @*/ 304af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 30529e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 30629e0b56fSBarry Smith { 30729e0b56fSBarry Smith int ierr; 30829e0b56fSBarry Smith Scalar mone = -1.0; 3098f99978dSLois Curfman McInnes SNES_LS *neP = (SNES_LS *) snes->data; 3108f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 31129e0b56fSBarry Smith 3123a40ed3dSBarry Smith PetscFunctionBegin; 313*8cbba510SBarry Smith *flag = 0; 31429e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 31529e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 3168f99978dSLois Curfman McInnes if (neP->CheckStep) { 3178f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 3188f99978dSLois Curfman McInnes } 31929e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 32029e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3213a40ed3dSBarry Smith PetscFunctionReturn(0); 32229e0b56fSBarry Smith } 32304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 32429e0b56fSBarry Smith #undef __FUNC__ 3255615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3264b828684SBarry Smith /*@C 327f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3285e76c431SBarry Smith 329c7afd0dbSLois Curfman McInnes Collective on SNES 330c7afd0dbSLois Curfman McInnes 3315e76c431SBarry Smith Input Parameters: 332c7afd0dbSLois Curfman McInnes + snes - nonlinear context 333af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3345e76c431SBarry Smith . x - current iterate 3355e76c431SBarry Smith . f - residual evaluated at x 3365e76c431SBarry Smith . y - search direction (contains new iterate on output) 3375e76c431SBarry Smith . w - work vector 338c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3395e76c431SBarry Smith 340393d2d9aSLois Curfman McInnes Output Parameters: 341c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3425e76c431SBarry Smith . y - new iterate (contains search direction on input) 3435e76c431SBarry Smith . gnorm - 2-norm of g 3445e76c431SBarry Smith . ynorm - 2-norm of search length 345c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 346fee21e36SBarry Smith 347c4a48953SLois Curfman McInnes Options Database Key: 348c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 349c4a48953SLois Curfman McInnes 3505e76c431SBarry Smith Notes: 3515e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3525e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3535e76c431SBarry Smith 35436851e7fSLois Curfman McInnes Level: advanced 35536851e7fSLois Curfman McInnes 35628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 35728ae5a21SLois Curfman McInnes 358af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3595e76c431SBarry Smith @*/ 360af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 361761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3625e76c431SBarry Smith { 363406556e6SLois Curfman McInnes /* 364406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 365406556e6SLois Curfman McInnes minimization problem: 366406556e6SLois Curfman McInnes min z(x): R^n -> R, 367406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 368406556e6SLois Curfman McInnes */ 369406556e6SLois Curfman McInnes 370ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 371ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 3723a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 3735e42470aSBarry Smith Scalar cinitslope, clambda; 3746b5873e3SBarry Smith #endif 3755e42470aSBarry Smith int ierr, count; 3765e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3775334005bSBarry Smith Scalar mone = -1.0, scale; 3788f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3795e76c431SBarry Smith 3803a40ed3dSBarry Smith PetscFunctionBegin; 3817857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 382761aaf1bSLois Curfman McInnes *flag = 0; 3835e76c431SBarry Smith alpha = neP->alpha; 3845e76c431SBarry Smith maxstep = neP->maxstep; 3855e76c431SBarry Smith steptol = neP->steptol; 3865e76c431SBarry Smith 387cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); 388a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 389a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 390a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 391a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 392a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 393ad922ac9SBarry Smith goto theend1; 39494a9d846SBarry Smith } 3955e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3965e42470aSBarry Smith scale = maxstep/(*ynorm); 3973a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 39892318cfeSSatish Balay PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale)); 3996b5873e3SBarry Smith #else 40094a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 4016b5873e3SBarry Smith #endif 402393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 4035e76c431SBarry Smith *ynorm = maxstep; 4045e76c431SBarry Smith } 4055e76c431SBarry Smith minlambda = steptol/(*ynorm); 406a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 4073a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 408a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 40992318cfeSSatish Balay initslope = PetscReal(cinitslope); 4105e42470aSBarry Smith #else 411a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 4125e42470aSBarry Smith #endif 4135e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4145e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4155e76c431SBarry Smith 416393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 4175334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 41878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 419cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); 420406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 421393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 42294a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 423ad922ac9SBarry Smith goto theend1; 4245e76c431SBarry Smith } 4255e76c431SBarry Smith 4265e76c431SBarry Smith /* Fit points with quadratic */ 4275e76c431SBarry Smith lambda = 1.0; count = 0; 428a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4295e76c431SBarry Smith lambdaprev = lambda; 4305e76c431SBarry Smith gnormprev = *gnorm; 431ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 432ddd12767SBarry Smith else lambda = lambdatemp; 433393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 434ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4353a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 436ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4375e42470aSBarry Smith #else 438ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4395e42470aSBarry Smith #endif 44078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 441cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 442406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 443393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 444ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 445ad922ac9SBarry Smith goto theend1; 4465e76c431SBarry Smith } 4475e76c431SBarry Smith 4485e76c431SBarry Smith /* Fit points with cubic */ 4495e76c431SBarry Smith count = 1; 4505e76c431SBarry Smith while (1) { 4515e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4522b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4532b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 454ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 455393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 456761aaf1bSLois Curfman McInnes *flag = -1; break; 4575e76c431SBarry Smith } 458406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 459406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 460ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4612b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4625e76c431SBarry Smith d = b*b - 3*a*initslope; 4635e76c431SBarry Smith if (d < 0.0) d = 0.0; 4645e76c431SBarry Smith if (a == 0.0) { 4655e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4665e76c431SBarry Smith } else { 4675e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4685e76c431SBarry Smith } 4695e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4705e76c431SBarry Smith lambdatemp = .5*lambda; 4715e76c431SBarry Smith } 4725e76c431SBarry Smith lambdaprev = lambda; 4735e76c431SBarry Smith gnormprev = *gnorm; 4745e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4755e76c431SBarry Smith lambda = .1*lambda; 4765e76c431SBarry Smith } 4775e76c431SBarry Smith else lambda = lambdatemp; 478393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 479ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 481ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 482393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4835e42470aSBarry Smith #else 484ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4855e42470aSBarry Smith #endif 48678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 487cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 488406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 489393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 490ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4912715565aSLois Curfman McInnes goto theend1; 4922b022350SLois Curfman McInnes } else { 4932b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 4945e76c431SBarry Smith } 4955e76c431SBarry Smith count++; 4965e76c431SBarry Smith } 497ad922ac9SBarry Smith theend1: 4988f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 4998f99978dSLois Curfman McInnes if (neP->CheckStep) { 5008f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 5018f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5028f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); 5038f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); 5048f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); 5058f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 5068f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr); 5078f99978dSLois Curfman McInnes } 5088f99978dSLois Curfman McInnes } 5097857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5103a40ed3dSBarry Smith PetscFunctionReturn(0); 5115e76c431SBarry Smith } 51204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5135615d1e5SSatish Balay #undef __FUNC__ 5145615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 5154b828684SBarry Smith /*@C 516f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5175e76c431SBarry Smith 518c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 519c7afd0dbSLois Curfman McInnes 5205e42470aSBarry Smith Input Parameters: 521c7afd0dbSLois Curfman McInnes + snes - the SNES context 522af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5235e42470aSBarry Smith . x - current iterate 5245e42470aSBarry Smith . f - residual evaluated at x 5255e42470aSBarry Smith . y - search direction (contains new iterate on output) 5265e42470aSBarry Smith . w - work vector 527c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5285e42470aSBarry Smith 529c4a48953SLois Curfman McInnes Output Parameters: 530c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5315e42470aSBarry Smith . y - new iterate (contains search direction on input) 5325e42470aSBarry Smith . gnorm - 2-norm of g 5335e42470aSBarry Smith . ynorm - 2-norm of search length 534c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 535fee21e36SBarry Smith 536c4a48953SLois Curfman McInnes Options Database Key: 537c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 538c4a48953SLois Curfman McInnes 5395e42470aSBarry Smith Notes: 540c7afd0dbSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method. 5415e42470aSBarry Smith 54236851e7fSLois Curfman McInnes Level: advanced 54336851e7fSLois Curfman McInnes 544f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 545f59ffdedSLois Curfman McInnes 546af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5475e42470aSBarry Smith @*/ 548af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 549761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5505e76c431SBarry Smith { 551406556e6SLois Curfman McInnes /* 552406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 553406556e6SLois Curfman McInnes minimization problem: 554406556e6SLois Curfman McInnes min z(x): R^n -> R, 555406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 556406556e6SLois Curfman McInnes */ 55740011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 5583a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 5595e42470aSBarry Smith Scalar cinitslope,clambda; 5606b5873e3SBarry Smith #endif 5615e42470aSBarry Smith int ierr, count; 5625e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 5635334005bSBarry Smith Scalar mone = -1.0,scale; 5648f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5655e76c431SBarry Smith 5663a40ed3dSBarry Smith PetscFunctionBegin; 5677857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 568761aaf1bSLois Curfman McInnes *flag = 0; 5695e42470aSBarry Smith alpha = neP->alpha; 5705e42470aSBarry Smith maxstep = neP->maxstep; 5715e42470aSBarry Smith steptol = neP->steptol; 5725e76c431SBarry Smith 573cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 574b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 57594a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 576b37302c6SLois Curfman McInnes *gnorm = fnorm; 577b37302c6SLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 578b37302c6SLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 579ad922ac9SBarry Smith goto theend2; 58094a9d846SBarry Smith } 5815e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5825e42470aSBarry Smith scale = maxstep/(*ynorm); 583393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 5845e42470aSBarry Smith *ynorm = maxstep; 5855e76c431SBarry Smith } 5865e42470aSBarry Smith minlambda = steptol/(*ynorm); 587a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 5883a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 589a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 59092318cfeSSatish Balay initslope = PetscReal(cinitslope); 5915e42470aSBarry Smith #else 592a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 5935e42470aSBarry Smith #endif 5945e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 5955e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 5965e42470aSBarry Smith 597393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 5985334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 59978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 600cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 601406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 602393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 60394a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 604ad922ac9SBarry Smith goto theend2; 6055e42470aSBarry Smith } 6065e42470aSBarry Smith 6075e42470aSBarry Smith /* Fit points with quadratic */ 6085e42470aSBarry Smith lambda = 1.0; count = 0; 6095e42470aSBarry Smith count = 1; 6105e42470aSBarry Smith while (1) { 6115e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 612981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 613981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 614ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 615393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 616761aaf1bSLois Curfman McInnes *flag = -1; break; 6175e42470aSBarry Smith } 618a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 6195e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 6205e42470aSBarry Smith lambda = .1*lambda; 6215e42470aSBarry Smith } else lambda = lambdatemp; 622393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 6235334005bSBarry Smith lambda = -lambda; 6243a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 625393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 6265e42470aSBarry Smith #else 627393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 6285e42470aSBarry Smith #endif 62978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 630cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 631406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 632393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 633981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 63406259719SBarry Smith break; 6355e42470aSBarry Smith } 6365e42470aSBarry Smith count++; 6375e42470aSBarry Smith } 638ad922ac9SBarry Smith theend2: 6398f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6408f99978dSLois Curfman McInnes if (neP->CheckStep) { 6418f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr); 6428f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6438f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); 6448f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr); 6458f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr); 6468f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr); 6478f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr); 6488f99978dSLois Curfman McInnes } 6498f99978dSLois Curfman McInnes } 6507857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6513a40ed3dSBarry Smith PetscFunctionReturn(0); 6525e76c431SBarry Smith } 65304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6545615d1e5SSatish Balay #undef __FUNC__ 655d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 656c9e6a524SLois Curfman McInnes /*@C 65777c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 658f63b844aSLois Curfman McInnes by the method SNES_EQ_LS. 6595e76c431SBarry Smith 6605e76c431SBarry Smith Input Parameters: 661c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 662af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 663c7afd0dbSLois Curfman McInnes - func - pointer to int function 6645e76c431SBarry Smith 665fee21e36SBarry Smith Collective on SNES 666fee21e36SBarry Smith 667c4a48953SLois Curfman McInnes Available Routines: 668c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 669f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 670f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 671af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6725e76c431SBarry Smith 673c4a48953SLois Curfman McInnes Options Database Keys: 674af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 675c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 676c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 677c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 678c4a48953SLois Curfman McInnes 6795e76c431SBarry Smith Calling sequence of func: 680bd208895SLois Curfman McInnes .vb 681af2d14d2SLois Curfman McInnes func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w, 682bd208895SLois Curfman McInnes double fnorm, double *ynorm, double *gnorm, *flag) 683bd208895SLois Curfman McInnes .ve 6845e76c431SBarry Smith 6855e76c431SBarry Smith Input parameters for func: 686c7afd0dbSLois Curfman McInnes + snes - nonlinear context 687af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 6885e76c431SBarry Smith . x - current iterate 6895e76c431SBarry Smith . f - residual evaluated at x 6905e76c431SBarry Smith . y - search direction (contains new iterate on output) 6915e76c431SBarry Smith . w - work vector 692c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6935e76c431SBarry Smith 6945e76c431SBarry Smith Output parameters for func: 695c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6965e76c431SBarry Smith . y - new iterate (contains search direction on input) 6975e76c431SBarry Smith . gnorm - 2-norm of g 6985e76c431SBarry Smith . ynorm - 2-norm of search length 699c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 700761aaf1bSLois Curfman McInnes on failure. 701f59ffdedSLois Curfman McInnes 70236851e7fSLois Curfman McInnes Level: advanced 70336851e7fSLois Curfman McInnes 704f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 705f59ffdedSLois Curfman McInnes 7068f99978dSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck() 7075e76c431SBarry Smith @*/ 708af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx) 7095e76c431SBarry Smith { 710af2d14d2SLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*); 71182bf6240SBarry Smith 7123a40ed3dSBarry Smith PetscFunctionBegin; 71394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 71482bf6240SBarry Smith if (f) { 715af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 71682bf6240SBarry Smith } 7173a40ed3dSBarry Smith PetscFunctionReturn(0); 7185e76c431SBarry Smith } 71904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 720fb2e594dSBarry Smith EXTERN_C_BEGIN 72182bf6240SBarry Smith #undef __FUNC__ 72282bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 723af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 724af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 72582bf6240SBarry Smith { 72682bf6240SBarry Smith PetscFunctionBegin; 72782bf6240SBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 728af2d14d2SLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 72982bf6240SBarry Smith PetscFunctionReturn(0); 73082bf6240SBarry Smith } 731fb2e594dSBarry Smith EXTERN_C_END 73204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 733c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 734c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck" 735c8dd0c0dSLois Curfman McInnes /*@C 7368f99978dSLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the new iterate computed 7378f99978dSLois Curfman McInnes by the line search routine in the Newton-based method SNES_EQ_LS. 738c8dd0c0dSLois Curfman McInnes 739c8dd0c0dSLois Curfman McInnes Input Parameters: 740c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 741c8dd0c0dSLois Curfman McInnes . func - pointer to int function 742c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 743c8dd0c0dSLois Curfman McInnes 744c8dd0c0dSLois Curfman McInnes Collective on SNES 745c8dd0c0dSLois Curfman McInnes 746c8dd0c0dSLois Curfman McInnes Calling sequence of func: 747c8dd0c0dSLois Curfman McInnes .vb 7488f99978dSLois Curfman McInnes func (SNES snes, void *lsctx, Vec x, PetscTruth *flag) 749c8dd0c0dSLois Curfman McInnes .ve 750c8dd0c0dSLois Curfman McInnes 751c8dd0c0dSLois Curfman McInnes Input parameters for func: 752c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 753c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7548f99978dSLois Curfman McInnes - x - current candidate iterate 755c8dd0c0dSLois Curfman McInnes 756c8dd0c0dSLois Curfman McInnes Output parameters for func: 757c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 758c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 759c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 760c8dd0c0dSLois Curfman McInnes 761c8dd0c0dSLois Curfman McInnes Level: advanced 762c8dd0c0dSLois Curfman McInnes 7638f99978dSLois Curfman McInnes Notes: 7648f99978dSLois Curfman McInnes The user-defined line search checking routine is available for 7658f99978dSLois Curfman McInnes use in conjunction with SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 7668f99978dSLois Curfman McInnes which (1) compute a candidate iterate u_{i+1}, (2) pass control to the 7678f99978dSLois Curfman McInnes checking routine, and then (3) compute the corresponding function f(u_{i+1}) 7688f99978dSLois Curfman McInnes with the (possibly altered) iterate u_{i+1}. 7698f99978dSLois Curfman McInnes 7708f99978dSLois Curfman McInnes The user-defined line search checking routine is also available for 7718f99978dSLois Curfman McInnes use in conjunction with SNESQuadraticLineSearch() and SNESCubicLineSearch(). 7728f99978dSLois Curfman McInnes These routines (1) compute a candidate iterate u_{i+1} as well as a 7738f99978dSLois Curfman McInnes candidate function f(u_{i+1}), (2) pass control to the checking routine, 7748f99978dSLois Curfman McInnes and then (3) force a re-evaluation of f(u_{i+1}) if any changes were made 7758f99978dSLois Curfman McInnes to the candidate iterate in the checking routine (as indicated by 7768f99978dSLois Curfman McInnes flag=PETSC_TRUE). The overhead of this re-evaluation can be costly, so 7778f99978dSLois Curfman McInnes this feature with caution. 7788f99978dSLois Curfman McInnes 779c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 780c8dd0c0dSLois Curfman McInnes 781c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 782c8dd0c0dSLois Curfman McInnes @*/ 7838f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 784c8dd0c0dSLois Curfman McInnes { 7858f99978dSLois Curfman McInnes int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 786c8dd0c0dSLois Curfman McInnes 787c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 788c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 789c8dd0c0dSLois Curfman McInnes if (f) { 790c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 791c8dd0c0dSLois Curfman McInnes } 792c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 793c8dd0c0dSLois Curfman McInnes } 794c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 795c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 796c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 797c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS" 7988f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 799c8dd0c0dSLois Curfman McInnes { 800c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 801c8dd0c0dSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 802c8dd0c0dSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 803c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 804c8dd0c0dSLois Curfman McInnes } 805c8dd0c0dSLois Curfman McInnes EXTERN_C_END 806c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 80704d965bbSLois Curfman McInnes /* 80804d965bbSLois Curfman McInnes SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 80982bf6240SBarry Smith 81004d965bbSLois Curfman McInnes Input Parameter: 81104d965bbSLois Curfman McInnes . snes - the SNES context 81204d965bbSLois Curfman McInnes 81304d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 81404d965bbSLois Curfman McInnes */ 8155615d1e5SSatish Balay #undef __FUNC__ 816d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 817f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8185e42470aSBarry Smith { 8195e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 8206daaf66cSBarry Smith 8213a40ed3dSBarry Smith PetscFunctionBegin; 82276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n"); 823af2d14d2SLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 82476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 82576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 82676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol); 8273a40ed3dSBarry Smith PetscFunctionReturn(0); 8285e42470aSBarry Smith } 82904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 83004d965bbSLois Curfman McInnes /* 831de49b36dSLois Curfman McInnes SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure. 83204d965bbSLois Curfman McInnes 83304d965bbSLois Curfman McInnes Input Parameters: 83404d965bbSLois Curfman McInnes . SNES - the SNES context 83504d965bbSLois Curfman McInnes . viewer - visualization context 83604d965bbSLois Curfman McInnes 83704d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 83804d965bbSLois Curfman McInnes */ 8395615d1e5SSatish Balay #undef __FUNC__ 840d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 841e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 842a935fc98SLois Curfman McInnes { 843a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 84419bcc07fSBarry Smith char *cstr; 84551695f54SLois Curfman McInnes int ierr; 84619bcc07fSBarry Smith ViewerType vtype; 847a935fc98SLois Curfman McInnes 8483a40ed3dSBarry Smith PetscFunctionBegin; 84919bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 8503f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 85119bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 85219bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 85319bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 85419bcc07fSBarry Smith else cstr = "unknown"; 8550ef38995SBarry Smith ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr); 8560ef38995SBarry Smith ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol); 8575cd90555SBarry Smith } else { 8585cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 85919bcc07fSBarry Smith } 8603a40ed3dSBarry Smith PetscFunctionReturn(0); 861a935fc98SLois Curfman McInnes } 86204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 86304d965bbSLois Curfman McInnes /* 86404d965bbSLois Curfman McInnes SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method. 86504d965bbSLois Curfman McInnes 86604d965bbSLois Curfman McInnes Input Parameter: 86704d965bbSLois Curfman McInnes . snes - the SNES context 86804d965bbSLois Curfman McInnes 86904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 87004d965bbSLois Curfman McInnes */ 8715615d1e5SSatish Balay #undef __FUNC__ 8725615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 873f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 8745e42470aSBarry Smith { 8755e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 8765e42470aSBarry Smith char ver[16]; 8775e42470aSBarry Smith double tmp; 87825cdf11fSBarry Smith int ierr,flg; 8795e42470aSBarry Smith 8803a40ed3dSBarry Smith PetscFunctionBegin; 88109d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr); 88225cdf11fSBarry Smith if (flg) { 8835e42470aSBarry Smith ls->alpha = tmp; 8845e42470aSBarry Smith } 88509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr); 88625cdf11fSBarry Smith if (flg) { 8875e42470aSBarry Smith ls->maxstep = tmp; 8885e42470aSBarry Smith } 88909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr); 89025cdf11fSBarry Smith if (flg) { 8915e42470aSBarry Smith ls->steptol = tmp; 8925e42470aSBarry Smith } 89309d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr); 89425cdf11fSBarry Smith if (flg) { 89548d91487SBarry Smith if (!PetscStrcmp(ver,"basic")) { 896af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 897b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"basicnonorms")) { 898af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 899b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"quadratic")) { 900af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 901b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"cubic")) { 902af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9035e42470aSBarry Smith } 904a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9055e42470aSBarry Smith } 9063a40ed3dSBarry Smith PetscFunctionReturn(0); 9075e42470aSBarry Smith } 90804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 90904d965bbSLois Curfman McInnes /* 91004d965bbSLois Curfman McInnes SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method, 91104d965bbSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 91204d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 91304d965bbSLois Curfman McInnes 91404d965bbSLois Curfman McInnes 91504d965bbSLois Curfman McInnes Input Parameter: 91604d965bbSLois Curfman McInnes . snes - the SNES context 91704d965bbSLois Curfman McInnes 91804d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 91904d965bbSLois Curfman McInnes */ 920fb2e594dSBarry Smith EXTERN_C_BEGIN 9215615d1e5SSatish Balay #undef __FUNC__ 9225615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 923f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9245e42470aSBarry Smith { 92582bf6240SBarry Smith int ierr; 9265e42470aSBarry Smith SNES_LS *neP; 9275e42470aSBarry Smith 9283a40ed3dSBarry Smith PetscFunctionBegin; 929a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 930a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 931a8c6a408SBarry Smith } 93282bf6240SBarry Smith 933f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 934f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 935f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 936f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 937f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 938f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 939f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9405baf8537SBarry Smith snes->nwork = 0; 9415e42470aSBarry Smith 9420452661fSBarry Smith neP = PetscNew(SNES_LS); CHKPTRQ(neP); 943ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 9445e42470aSBarry Smith snes->data = (void *) neP; 9455e42470aSBarry Smith neP->alpha = 1.e-4; 9465e42470aSBarry Smith neP->maxstep = 1.e8; 9475e42470aSBarry Smith neP->steptol = 1.e-12; 9485e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 949c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 950c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 951c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 95282bf6240SBarry Smith 95394d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 954e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 955c8dd0c0dSLois Curfman McInnes ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 956c8dd0c0dSLois Curfman McInnes (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 95782bf6240SBarry Smith 9583a40ed3dSBarry Smith PetscFunctionReturn(0); 9595e42470aSBarry Smith } 960fb2e594dSBarry Smith EXTERN_C_END 96182bf6240SBarry Smith 96282bf6240SBarry Smith 96382bf6240SBarry Smith 964