1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*84c569c9SBarry Smith static char vcid[] = "$Id: ls.c,v 1.123 1999/02/09 23:32:19 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; 69*84c569c9SBarry Smith int maxits, i, ierr, lits, lsfail; 70112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 71*84c569c9SBarry 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); 85*84c569c9SBarry Smith snes->iter = 0; 865e42470aSBarry Smith snes->norm = fnorm; 8706d0569aSBarry Smith PetscAMSGrantAccess(snes); 88*84c569c9SBarry 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); 110ddd12767SBarry Smith ierr = (*neP->LineSearch)(snes,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); 119*84c569c9SBarry Smith snes->iter = i+1; 1205e42470aSBarry Smith snes->norm = fnorm; 12106d0569aSBarry Smith PetscAMSGrantAccess(snes); 122*84c569c9SBarry 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 2105e76c431SBarry Smith . x - current iterate 2115e76c431SBarry Smith . f - residual evaluated at x 2125e76c431SBarry Smith . y - search direction (contains new iterate on output) 2135e76c431SBarry Smith . w - work vector 214c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2155e76c431SBarry Smith 216c4a48953SLois Curfman McInnes Output Parameters: 217c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2185e76c431SBarry Smith . y - new iterate (contains search direction on input) 2195e76c431SBarry Smith . gnorm - 2-norm of g 2205e76c431SBarry Smith . ynorm - 2-norm of search length 221c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 222fee21e36SBarry Smith 223c4a48953SLois Curfman McInnes Options Database Key: 224c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 225c4a48953SLois Curfman McInnes 22636851e7fSLois Curfman McInnes Level: advanced 22736851e7fSLois Curfman McInnes 22828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 22928ae5a21SLois Curfman McInnes 230f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 23177c4ece6SBarry Smith SNESSetLineSearch() 2325e76c431SBarry Smith @*/ 2335e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 234761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 2355e76c431SBarry Smith { 2365e42470aSBarry Smith int ierr; 2375334005bSBarry Smith Scalar mone = -1.0; 2385334005bSBarry Smith 2393a40ed3dSBarry Smith PetscFunctionBegin; 240761aaf1bSLois Curfman McInnes *flag = 0; 2417857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 242cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 243ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 244ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 245cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 2467857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2473a40ed3dSBarry Smith PetscFunctionReturn(0); 2485e76c431SBarry Smith } 24904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2505615d1e5SSatish Balay #undef __FUNC__ 25129e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 25282bf6240SBarry Smith 25329e0b56fSBarry Smith /*@C 25429e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 25529e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 25629e0b56fSBarry Smith even compute the norm of the function or search direction; this 25729e0b56fSBarry Smith is intended only when you know the full step is fine and are 25829e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 259c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 260c7afd0dbSLois Curfman McInnes 261c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 26229e0b56fSBarry Smith 26329e0b56fSBarry Smith Input Parameters: 264c7afd0dbSLois Curfman McInnes + snes - nonlinear context 26529e0b56fSBarry Smith . x - current iterate 26629e0b56fSBarry Smith . f - residual evaluated at x 26729e0b56fSBarry Smith . y - search direction (contains new iterate on output) 26829e0b56fSBarry Smith . w - work vector 269c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 27029e0b56fSBarry Smith 27129e0b56fSBarry Smith Output Parameters: 272c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 27329e0b56fSBarry Smith . gnorm - not changed 27429e0b56fSBarry Smith . ynorm - not changed 275c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 276fee21e36SBarry Smith 27729e0b56fSBarry Smith Options Database Key: 278c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 27929e0b56fSBarry Smith 28036851e7fSLois Curfman McInnes Level: advanced 28136851e7fSLois Curfman McInnes 28229e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 28329e0b56fSBarry Smith 28429e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 28529e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 28629e0b56fSBarry Smith @*/ 28729e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 28829e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 28929e0b56fSBarry Smith { 29029e0b56fSBarry Smith int ierr; 29129e0b56fSBarry Smith Scalar mone = -1.0; 29229e0b56fSBarry Smith 2933a40ed3dSBarry Smith PetscFunctionBegin; 29429e0b56fSBarry Smith *flag = 0; 29529e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 29629e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 29729e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 29829e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2993a40ed3dSBarry Smith PetscFunctionReturn(0); 30029e0b56fSBarry Smith } 30104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 30229e0b56fSBarry Smith #undef __FUNC__ 3035615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3044b828684SBarry Smith /*@C 305f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3065e76c431SBarry Smith 307c7afd0dbSLois Curfman McInnes Collective on SNES 308c7afd0dbSLois Curfman McInnes 3095e76c431SBarry Smith Input Parameters: 310c7afd0dbSLois Curfman McInnes + snes - nonlinear context 3115e76c431SBarry Smith . x - current iterate 3125e76c431SBarry Smith . f - residual evaluated at x 3135e76c431SBarry Smith . y - search direction (contains new iterate on output) 3145e76c431SBarry Smith . w - work vector 315c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3165e76c431SBarry Smith 317393d2d9aSLois Curfman McInnes Output Parameters: 318c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3195e76c431SBarry Smith . y - new iterate (contains search direction on input) 3205e76c431SBarry Smith . gnorm - 2-norm of g 3215e76c431SBarry Smith . ynorm - 2-norm of search length 322c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 323fee21e36SBarry Smith 324c4a48953SLois Curfman McInnes Options Database Key: 325c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 326c4a48953SLois Curfman McInnes 3275e76c431SBarry Smith Notes: 3285e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3295e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3305e76c431SBarry Smith 33136851e7fSLois Curfman McInnes Level: advanced 33236851e7fSLois Curfman McInnes 33328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 33428ae5a21SLois Curfman McInnes 33577c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 3365e76c431SBarry Smith @*/ 3375e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w, 338761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3395e76c431SBarry Smith { 340406556e6SLois Curfman McInnes /* 341406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 342406556e6SLois Curfman McInnes minimization problem: 343406556e6SLois Curfman McInnes min z(x): R^n -> R, 344406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 345406556e6SLois Curfman McInnes */ 346406556e6SLois Curfman McInnes 347ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 348ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 3493a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 3505e42470aSBarry Smith Scalar cinitslope, clambda; 3516b5873e3SBarry Smith #endif 3525e42470aSBarry Smith int ierr, count; 3535e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3545334005bSBarry Smith Scalar mone = -1.0,scale; 3555e76c431SBarry Smith 3563a40ed3dSBarry Smith PetscFunctionBegin; 3577857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 358761aaf1bSLois Curfman McInnes *flag = 0; 3595e76c431SBarry Smith alpha = neP->alpha; 3605e76c431SBarry Smith maxstep = neP->maxstep; 3615e76c431SBarry Smith steptol = neP->steptol; 3625e76c431SBarry Smith 363cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); 364a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 365a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 366a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 367a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 368a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 369ad922ac9SBarry Smith goto theend1; 37094a9d846SBarry Smith } 3715e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3725e42470aSBarry Smith scale = maxstep/(*ynorm); 3733a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 37492318cfeSSatish Balay PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale)); 3756b5873e3SBarry Smith #else 37694a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 3776b5873e3SBarry Smith #endif 378393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 3795e76c431SBarry Smith *ynorm = maxstep; 3805e76c431SBarry Smith } 3815e76c431SBarry Smith minlambda = steptol/(*ynorm); 382a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 3833a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 384a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 38592318cfeSSatish Balay initslope = PetscReal(cinitslope); 3865e42470aSBarry Smith #else 387a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 3885e42470aSBarry Smith #endif 3895e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 3905e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 3915e76c431SBarry Smith 392393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 3935334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 39478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 395cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); 396406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 397393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 39894a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 399ad922ac9SBarry Smith goto theend1; 4005e76c431SBarry Smith } 4015e76c431SBarry Smith 4025e76c431SBarry Smith /* Fit points with quadratic */ 4035e76c431SBarry Smith lambda = 1.0; count = 0; 404a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4055e76c431SBarry Smith lambdaprev = lambda; 4065e76c431SBarry Smith gnormprev = *gnorm; 407ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 408ddd12767SBarry Smith else lambda = lambdatemp; 409393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 410ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4113a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 412ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4135e42470aSBarry Smith #else 414ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4155e42470aSBarry Smith #endif 41678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 417cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 418406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 419393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 420ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 421ad922ac9SBarry Smith goto theend1; 4225e76c431SBarry Smith } 4235e76c431SBarry Smith 4245e76c431SBarry Smith /* Fit points with cubic */ 4255e76c431SBarry Smith count = 1; 4265e76c431SBarry Smith while (1) { 4275e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4282b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4292b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 430ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 431393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 432761aaf1bSLois Curfman McInnes *flag = -1; break; 4335e76c431SBarry Smith } 434406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 435406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 436ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4372b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4385e76c431SBarry Smith d = b*b - 3*a*initslope; 4395e76c431SBarry Smith if (d < 0.0) d = 0.0; 4405e76c431SBarry Smith if (a == 0.0) { 4415e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4425e76c431SBarry Smith } else { 4435e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4445e76c431SBarry Smith } 4455e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4465e76c431SBarry Smith lambdatemp = .5*lambda; 4475e76c431SBarry Smith } 4485e76c431SBarry Smith lambdaprev = lambda; 4495e76c431SBarry Smith gnormprev = *gnorm; 4505e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4515e76c431SBarry Smith lambda = .1*lambda; 4525e76c431SBarry Smith } 4535e76c431SBarry Smith else lambda = lambdatemp; 454393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 455ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4563a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 457ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 458393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4595e42470aSBarry Smith #else 460ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4615e42470aSBarry Smith #endif 46278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 463cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 464406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 465393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 466ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4672715565aSLois Curfman McInnes goto theend1; 4682b022350SLois Curfman McInnes } else { 4692b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 4705e76c431SBarry Smith } 4715e76c431SBarry Smith count++; 4725e76c431SBarry Smith } 473ad922ac9SBarry Smith theend1: 4747857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4753a40ed3dSBarry Smith PetscFunctionReturn(0); 4765e76c431SBarry Smith } 47704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4785615d1e5SSatish Balay #undef __FUNC__ 4795615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 4804b828684SBarry Smith /*@C 481f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 4825e76c431SBarry Smith 483c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 484c7afd0dbSLois Curfman McInnes 4855e42470aSBarry Smith Input Parameters: 486c7afd0dbSLois Curfman McInnes + snes - the SNES context 4875e42470aSBarry Smith . x - current iterate 4885e42470aSBarry Smith . f - residual evaluated at x 4895e42470aSBarry Smith . y - search direction (contains new iterate on output) 4905e42470aSBarry Smith . w - work vector 491c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4925e42470aSBarry Smith 493c4a48953SLois Curfman McInnes Output Parameters: 494c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4955e42470aSBarry Smith . y - new iterate (contains search direction on input) 4965e42470aSBarry Smith . gnorm - 2-norm of g 4975e42470aSBarry Smith . ynorm - 2-norm of search length 498c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 499fee21e36SBarry Smith 500c4a48953SLois Curfman McInnes Options Database Key: 501c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 502c4a48953SLois Curfman McInnes 5035e42470aSBarry Smith Notes: 504c7afd0dbSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method. 5055e42470aSBarry Smith 50636851e7fSLois Curfman McInnes Level: advanced 50736851e7fSLois Curfman McInnes 508f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 509f59ffdedSLois Curfman McInnes 51077c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 5115e42470aSBarry Smith @*/ 5125e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 513761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5145e76c431SBarry Smith { 515406556e6SLois Curfman McInnes /* 516406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 517406556e6SLois Curfman McInnes minimization problem: 518406556e6SLois Curfman McInnes min z(x): R^n -> R, 519406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 520406556e6SLois Curfman McInnes */ 52140011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 5223a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 5235e42470aSBarry Smith Scalar cinitslope,clambda; 5246b5873e3SBarry Smith #endif 5255e42470aSBarry Smith int ierr,count; 5265e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 5275334005bSBarry Smith Scalar mone = -1.0,scale; 5285e76c431SBarry Smith 5293a40ed3dSBarry Smith PetscFunctionBegin; 5307857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 531761aaf1bSLois Curfman McInnes *flag = 0; 5325e42470aSBarry Smith alpha = neP->alpha; 5335e42470aSBarry Smith maxstep = neP->maxstep; 5345e42470aSBarry Smith steptol = neP->steptol; 5355e76c431SBarry Smith 536cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 537b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 53894a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 539b37302c6SLois Curfman McInnes *gnorm = fnorm; 540b37302c6SLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 541b37302c6SLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 542ad922ac9SBarry Smith goto theend2; 54394a9d846SBarry Smith } 5445e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5455e42470aSBarry Smith scale = maxstep/(*ynorm); 546393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 5475e42470aSBarry Smith *ynorm = maxstep; 5485e76c431SBarry Smith } 5495e42470aSBarry Smith minlambda = steptol/(*ynorm); 550a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 5513a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 552a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 55392318cfeSSatish Balay initslope = PetscReal(cinitslope); 5545e42470aSBarry Smith #else 555a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 5565e42470aSBarry Smith #endif 5575e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 5585e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 5595e42470aSBarry Smith 560393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 5615334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 56278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 563cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 564406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 565393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 56694a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 567ad922ac9SBarry Smith goto theend2; 5685e42470aSBarry Smith } 5695e42470aSBarry Smith 5705e42470aSBarry Smith /* Fit points with quadratic */ 5715e42470aSBarry Smith lambda = 1.0; count = 0; 5725e42470aSBarry Smith count = 1; 5735e42470aSBarry Smith while (1) { 5745e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 575981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 576981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 577ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 578393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 579761aaf1bSLois Curfman McInnes *flag = -1; break; 5805e42470aSBarry Smith } 581a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5825e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 5835e42470aSBarry Smith lambda = .1*lambda; 5845e42470aSBarry Smith } else lambda = lambdatemp; 585393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 5865334005bSBarry Smith lambda = -lambda; 5873a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 588393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 5895e42470aSBarry Smith #else 590393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 5915e42470aSBarry Smith #endif 59278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 593cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 594406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 595393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 596981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 59706259719SBarry Smith break; 5985e42470aSBarry Smith } 5995e42470aSBarry Smith count++; 6005e42470aSBarry Smith } 601ad922ac9SBarry Smith theend2: 6027857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6033a40ed3dSBarry Smith PetscFunctionReturn(0); 6045e76c431SBarry Smith } 60504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6065615d1e5SSatish Balay #undef __FUNC__ 607d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 608c9e6a524SLois Curfman McInnes /*@C 60977c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 610f63b844aSLois Curfman McInnes by the method SNES_EQ_LS. 6115e76c431SBarry Smith 6125e76c431SBarry Smith Input Parameters: 613c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 614c7afd0dbSLois Curfman McInnes - func - pointer to int function 6155e76c431SBarry Smith 616fee21e36SBarry Smith Collective on SNES 617fee21e36SBarry Smith 618c4a48953SLois Curfman McInnes Available Routines: 619c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 620f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 621f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 622c7afd0dbSLois Curfman McInnes - SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel) 6235e76c431SBarry Smith 624c4a48953SLois Curfman McInnes Options Database Keys: 625c7afd0dbSLois Curfman McInnes + -snes_eq_ls [basic,quadratic,cubic] - Selects line search 626c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 627c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 628c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 629c4a48953SLois Curfman McInnes 6305e76c431SBarry Smith Calling sequence of func: 631bd208895SLois Curfman McInnes .vb 632bd208895SLois Curfman McInnes func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 633bd208895SLois Curfman McInnes double fnorm, double *ynorm, double *gnorm, *flag) 634bd208895SLois Curfman McInnes .ve 6355e76c431SBarry Smith 6365e76c431SBarry Smith Input parameters for func: 637c7afd0dbSLois Curfman McInnes + snes - nonlinear context 6385e76c431SBarry Smith . x - current iterate 6395e76c431SBarry Smith . f - residual evaluated at x 6405e76c431SBarry Smith . y - search direction (contains new iterate on output) 6415e76c431SBarry Smith . w - work vector 642c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6435e76c431SBarry Smith 6445e76c431SBarry Smith Output parameters for func: 645c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6465e76c431SBarry Smith . y - new iterate (contains search direction on input) 6475e76c431SBarry Smith . gnorm - 2-norm of g 6485e76c431SBarry Smith . ynorm - 2-norm of search length 649c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 650761aaf1bSLois Curfman McInnes on failure. 651f59ffdedSLois Curfman McInnes 65236851e7fSLois Curfman McInnes Level: advanced 65336851e7fSLois Curfman McInnes 654f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 655f59ffdedSLois Curfman McInnes 656f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch() 6575e76c431SBarry Smith @*/ 658*84c569c9SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*)) 6595e76c431SBarry Smith { 660415aef20SBarry Smith int ierr, (*f)(SNES,int (*)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*)); 66182bf6240SBarry Smith 6623a40ed3dSBarry Smith PetscFunctionBegin; 66394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 66482bf6240SBarry Smith if (f) { 66582bf6240SBarry Smith ierr = (*f)(snes,func);CHKERRQ(ierr); 66682bf6240SBarry Smith } 6673a40ed3dSBarry Smith PetscFunctionReturn(0); 6685e76c431SBarry Smith } 66904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 670fb2e594dSBarry Smith EXTERN_C_BEGIN 67182bf6240SBarry Smith #undef __FUNC__ 67282bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 67382bf6240SBarry Smith int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 67482bf6240SBarry Smith double,double*,double*,int*)) 67582bf6240SBarry Smith { 67682bf6240SBarry Smith PetscFunctionBegin; 67782bf6240SBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 67882bf6240SBarry Smith PetscFunctionReturn(0); 67982bf6240SBarry Smith } 680fb2e594dSBarry Smith EXTERN_C_END 68104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 68204d965bbSLois Curfman McInnes /* 68304d965bbSLois Curfman McInnes SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 68482bf6240SBarry Smith 68504d965bbSLois Curfman McInnes Input Parameter: 68604d965bbSLois Curfman McInnes . snes - the SNES context 68704d965bbSLois Curfman McInnes 68804d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 68904d965bbSLois Curfman McInnes */ 6905615d1e5SSatish Balay #undef __FUNC__ 691d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 692f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 6935e42470aSBarry Smith { 6945e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 6956daaf66cSBarry Smith 6963a40ed3dSBarry Smith PetscFunctionBegin; 69776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n"); 69876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [basic,quadratic,cubic]\n",p); 69976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 70076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 70176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol); 7023a40ed3dSBarry Smith PetscFunctionReturn(0); 7035e42470aSBarry Smith } 70404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 70504d965bbSLois Curfman McInnes /* 706de49b36dSLois Curfman McInnes SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure. 70704d965bbSLois Curfman McInnes 70804d965bbSLois Curfman McInnes Input Parameters: 70904d965bbSLois Curfman McInnes . SNES - the SNES context 71004d965bbSLois Curfman McInnes . viewer - visualization context 71104d965bbSLois Curfman McInnes 71204d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 71304d965bbSLois Curfman McInnes */ 7145615d1e5SSatish Balay #undef __FUNC__ 715d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 716e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 717a935fc98SLois Curfman McInnes { 718a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 71919bcc07fSBarry Smith char *cstr; 72051695f54SLois Curfman McInnes int ierr; 72119bcc07fSBarry Smith ViewerType vtype; 722a935fc98SLois Curfman McInnes 7233a40ed3dSBarry Smith PetscFunctionBegin; 72419bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 7253f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 72619bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 72719bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 72819bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 72919bcc07fSBarry Smith else cstr = "unknown"; 7300ef38995SBarry Smith ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr); 7310ef38995SBarry Smith ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol); 7325cd90555SBarry Smith } else { 7335cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 73419bcc07fSBarry Smith } 7353a40ed3dSBarry Smith PetscFunctionReturn(0); 736a935fc98SLois Curfman McInnes } 73704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 73804d965bbSLois Curfman McInnes /* 73904d965bbSLois Curfman McInnes SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method. 74004d965bbSLois Curfman McInnes 74104d965bbSLois Curfman McInnes Input Parameter: 74204d965bbSLois Curfman McInnes . snes - the SNES context 74304d965bbSLois Curfman McInnes 74404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 74504d965bbSLois Curfman McInnes */ 7465615d1e5SSatish Balay #undef __FUNC__ 7475615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 748f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 7495e42470aSBarry Smith { 7505e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 7515e42470aSBarry Smith char ver[16]; 7525e42470aSBarry Smith double tmp; 75325cdf11fSBarry Smith int ierr,flg; 7545e42470aSBarry Smith 7553a40ed3dSBarry Smith PetscFunctionBegin; 75609d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr); 75725cdf11fSBarry Smith if (flg) { 7585e42470aSBarry Smith ls->alpha = tmp; 7595e42470aSBarry Smith } 76009d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr); 76125cdf11fSBarry Smith if (flg) { 7625e42470aSBarry Smith ls->maxstep = tmp; 7635e42470aSBarry Smith } 76409d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr); 76525cdf11fSBarry Smith if (flg) { 7665e42470aSBarry Smith ls->steptol = tmp; 7675e42470aSBarry Smith } 76809d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr); 76925cdf11fSBarry Smith if (flg) { 77048d91487SBarry Smith if (!PetscStrcmp(ver,"basic")) { 771b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESNoLineSearch);CHKERRQ(ierr); 772b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"basicnonorms")) { 773b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);CHKERRQ(ierr); 774b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"quadratic")) { 775b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch);CHKERRQ(ierr); 776b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"cubic")) { 777b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESCubicLineSearch);CHKERRQ(ierr); 7785e42470aSBarry Smith } 779a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 7805e42470aSBarry Smith } 7813a40ed3dSBarry Smith PetscFunctionReturn(0); 7825e42470aSBarry Smith } 78304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 78404d965bbSLois Curfman McInnes /* 78504d965bbSLois Curfman McInnes SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method, 78604d965bbSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 78704d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 78804d965bbSLois Curfman McInnes 78904d965bbSLois Curfman McInnes 79004d965bbSLois Curfman McInnes Input Parameter: 79104d965bbSLois Curfman McInnes . snes - the SNES context 79204d965bbSLois Curfman McInnes 79304d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 79404d965bbSLois Curfman McInnes */ 795fb2e594dSBarry Smith EXTERN_C_BEGIN 7965615d1e5SSatish Balay #undef __FUNC__ 7975615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 798f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 7995e42470aSBarry Smith { 80082bf6240SBarry Smith int ierr; 8015e42470aSBarry Smith SNES_LS *neP; 8025e42470aSBarry Smith 8033a40ed3dSBarry Smith PetscFunctionBegin; 804a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 805a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 806a8c6a408SBarry Smith } 80782bf6240SBarry Smith 808f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 809f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 810f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 811f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 812f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 813f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 814f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 8155baf8537SBarry Smith snes->nwork = 0; 8165e42470aSBarry Smith 8170452661fSBarry Smith neP = PetscNew(SNES_LS); CHKPTRQ(neP); 818ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 8195e42470aSBarry Smith snes->data = (void *) neP; 8205e42470aSBarry Smith neP->alpha = 1.e-4; 8215e42470aSBarry Smith neP->maxstep = 1.e8; 8225e42470aSBarry Smith neP->steptol = 1.e-12; 8235e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 82482bf6240SBarry Smith 82594d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 826e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 82782bf6240SBarry Smith 8283a40ed3dSBarry Smith PetscFunctionReturn(0); 8295e42470aSBarry Smith } 830fb2e594dSBarry Smith EXTERN_C_END 83182bf6240SBarry Smith 83282bf6240SBarry Smith 83382bf6240SBarry Smith 834