1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*c7afd0dbSLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.108 1998/04/22 23:59:55 curfman Exp curfman $"; 35e76c431SBarry Smith #endif 45e76c431SBarry Smith 55e76c431SBarry Smith #include <math.h> 670f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 7b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 85e42470aSBarry Smith 904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 1004d965bbSLois Curfman McInnes 1104d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 1204d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 1304d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 1404d965bbSLois Curfman McInnes respectively. 1504d965bbSLois Curfman McInnes 16fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 1704d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 1804d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 19fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 2004d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 2104d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 2204d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 2304d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 2404d965bbSLois Curfman McInnes These routines are actually called via the common user interface 2504d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 2604d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 2704d965bbSLois Curfman McInnes for all nonlinear solvers. 2804d965bbSLois Curfman McInnes 2904d965bbSLois Curfman McInnes Another key routine is: 3004d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 3104d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 3204d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 3304d965bbSLois Curfman McInnes SNESSolve() if necessary. 3404d965bbSLois Curfman McInnes 3504d965bbSLois Curfman McInnes Additional basic routines are: 3604d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 3704d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 3804d965bbSLois Curfman McInnes have actually been used. 3904d965bbSLois Curfman McInnes These are called by application codes via the interface routines 4004d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 4104d965bbSLois Curfman McInnes 4204d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 4304d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 4404d965bbSLois Curfman McInnes above description applies to these categories also. 4504d965bbSLois Curfman McInnes 4604d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 475e42470aSBarry Smith /* 4804d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 4904d965bbSLois Curfman McInnes method with a line search. 505e76c431SBarry Smith 5104d965bbSLois Curfman McInnes Input Parameters: 5204d965bbSLois Curfman McInnes . snes - the SNES context 535e76c431SBarry Smith 5404d965bbSLois Curfman McInnes Output Parameter: 55c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 5604d965bbSLois Curfman McInnes 5704d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 585e76c431SBarry Smith 595e76c431SBarry Smith Notes: 605e76c431SBarry Smith This implements essentially a truncated Newton method with a 615e76c431SBarry Smith line search. By default a cubic backtracking line search 625e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 635e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 64393d2d9aSLois Curfman McInnes and Schnabel. 655e76c431SBarry Smith */ 665615d1e5SSatish Balay #undef __FUNC__ 675615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS" 68f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 695e76c431SBarry Smith { 705e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 71761aaf1bSLois Curfman McInnes int maxits, i, history_len, ierr, lits, lsfail; 72112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 736b5873e3SBarry Smith double fnorm, gnorm, xnorm, ynorm, *history; 745e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 755e76c431SBarry Smith 763a40ed3dSBarry Smith PetscFunctionBegin; 775e42470aSBarry Smith history = snes->conv_hist; /* convergence history */ 7851979daaSLois Curfman McInnes history_len = snes->conv_hist_size; /* convergence history length */ 795e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 805e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 8139e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 825e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 835e42470aSBarry Smith G = snes->work[1]; 845e42470aSBarry Smith W = snes->work[2]; 855e76c431SBarry Smith 8625ed9cd1SBarry Smith snes->iter = 0; 875334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr); /* F(X) */ 88cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr); /* fnorm <- ||F|| */ 895e42470aSBarry Smith snes->norm = fnorm; 9051979daaSLois Curfman McInnes if (history) history[0] = fnorm; 9194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 925e76c431SBarry Smith 933a40ed3dSBarry Smith if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);} 9494a9d846SBarry Smith 95d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 96d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 97d034289bSBarry Smith 985e76c431SBarry Smith for ( i=0; i<maxits; i++ ) { 995e42470aSBarry Smith snes->iter = i+1; 1005e76c431SBarry Smith 101ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1025334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr); 1035334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr); 10478b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr); 1057a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(lits); 106af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 107ea4d3ed3SLois Curfman McInnes 108ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 109ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 110ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 111ea4d3ed3SLois Curfman McInnes */ 11281b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr); 113ddd12767SBarry Smith ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr); 114af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 115761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 1165e76c431SBarry Smith 11739e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 11839e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1195e76c431SBarry Smith fnorm = gnorm; 1205e76c431SBarry Smith 1215e42470aSBarry Smith snes->norm = fnorm; 1225e76c431SBarry Smith if (history && history_len > i+1) history[i+1] = fnorm; 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 } 14251979daaSLois Curfman McInnes if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1; 1435e42470aSBarry Smith *outits = i+1; 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 1455e76c431SBarry Smith } 14604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 14704d965bbSLois Curfman McInnes /* 14804d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 14904d965bbSLois Curfman McInnes of the SNES_EQ_LS nonlinear solver. 15004d965bbSLois Curfman McInnes 15104d965bbSLois Curfman McInnes Input Parameter: 15204d965bbSLois Curfman McInnes . snes - the SNES context 15304d965bbSLois Curfman McInnes . x - the solution vector 15404d965bbSLois Curfman McInnes 15504d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 15604d965bbSLois Curfman McInnes 15704d965bbSLois Curfman McInnes Notes: 15804d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 15904d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 16004d965bbSLois Curfman McInnes the call to SNESSolve(). 16104d965bbSLois Curfman McInnes */ 1625615d1e5SSatish Balay #undef __FUNC__ 1635615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS" 164f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1655e76c431SBarry Smith { 1665e42470aSBarry Smith int ierr; 1673a40ed3dSBarry Smith 1683a40ed3dSBarry Smith PetscFunctionBegin; 16981b6cf68SLois Curfman McInnes snes->nwork = 4; 170d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1715e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 17281b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1733a40ed3dSBarry Smith PetscFunctionReturn(0); 1745e76c431SBarry Smith } 17504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 17604d965bbSLois Curfman McInnes /* 17704d965bbSLois Curfman McInnes SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created 17804d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 17904d965bbSLois Curfman McInnes 18004d965bbSLois Curfman McInnes Input Parameter: 18104d965bbSLois Curfman McInnes . snes - the SNES context 18204d965bbSLois Curfman McInnes 18304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetDestroy() 18404d965bbSLois Curfman McInnes */ 1855615d1e5SSatish Balay #undef __FUNC__ 186d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS" 187e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1885e76c431SBarry Smith { 189393d2d9aSLois Curfman McInnes int ierr; 1903a40ed3dSBarry Smith 1913a40ed3dSBarry Smith PetscFunctionBegin; 1925baf8537SBarry Smith if (snes->nwork) { 1934b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr); 19421c89e3eSBarry Smith } 1950452661fSBarry Smith PetscFree(snes->data); 1963a40ed3dSBarry Smith PetscFunctionReturn(0); 1975e76c431SBarry Smith } 19804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1995615d1e5SSatish Balay #undef __FUNC__ 2005615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch" 20182bf6240SBarry Smith 2024b828684SBarry Smith /*@C 2035e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2045e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2055e76c431SBarry Smith to serve as a template and is not recommended for general use. 2065e76c431SBarry Smith 207*c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 208*c7afd0dbSLois Curfman McInnes 2095e76c431SBarry Smith Input Parameters: 210*c7afd0dbSLois Curfman McInnes + snes - nonlinear context 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 215*c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2165e76c431SBarry Smith 217c4a48953SLois Curfman McInnes Output Parameters: 218*c7afd0dbSLois 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 222*c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 223fee21e36SBarry Smith 224c4a48953SLois Curfman McInnes Options Database Key: 225*c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 226c4a48953SLois Curfman McInnes 22728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 22828ae5a21SLois Curfman McInnes 229f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 23077c4ece6SBarry Smith SNESSetLineSearch() 2315e76c431SBarry Smith @*/ 2325e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 233761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 2345e76c431SBarry Smith { 2355e42470aSBarry Smith int ierr; 2365334005bSBarry Smith Scalar mone = -1.0; 2375334005bSBarry Smith 2383a40ed3dSBarry Smith PetscFunctionBegin; 239761aaf1bSLois Curfman McInnes *flag = 0; 2407857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 241cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 242ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 243ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 244cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 2457857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2475e76c431SBarry Smith } 24804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2495615d1e5SSatish Balay #undef __FUNC__ 25029e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 25182bf6240SBarry Smith 25229e0b56fSBarry Smith /*@C 25329e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 25429e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 25529e0b56fSBarry Smith even compute the norm of the function or search direction; this 25629e0b56fSBarry Smith is intended only when you know the full step is fine and are 25729e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 258*c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 259*c7afd0dbSLois Curfman McInnes 260*c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 26129e0b56fSBarry Smith 26229e0b56fSBarry Smith Input Parameters: 263*c7afd0dbSLois Curfman McInnes + snes - nonlinear context 26429e0b56fSBarry Smith . x - current iterate 26529e0b56fSBarry Smith . f - residual evaluated at x 26629e0b56fSBarry Smith . y - search direction (contains new iterate on output) 26729e0b56fSBarry Smith . w - work vector 268*c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 26929e0b56fSBarry Smith 27029e0b56fSBarry Smith Output Parameters: 271*c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 27229e0b56fSBarry Smith . gnorm - not changed 27329e0b56fSBarry Smith . ynorm - not changed 274*c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 275fee21e36SBarry Smith 27629e0b56fSBarry Smith Options Database Key: 277*c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 27829e0b56fSBarry Smith 27929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 28029e0b56fSBarry Smith 28129e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 28229e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 28329e0b56fSBarry Smith @*/ 28429e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 28529e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 28629e0b56fSBarry Smith { 28729e0b56fSBarry Smith int ierr; 28829e0b56fSBarry Smith Scalar mone = -1.0; 28929e0b56fSBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 29129e0b56fSBarry Smith *flag = 0; 29229e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 29329e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 29429e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 29529e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2963a40ed3dSBarry Smith PetscFunctionReturn(0); 29729e0b56fSBarry Smith } 29804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 29929e0b56fSBarry Smith #undef __FUNC__ 3005615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 3014b828684SBarry Smith /*@C 302f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3035e76c431SBarry Smith 304*c7afd0dbSLois Curfman McInnes Collective on SNES 305*c7afd0dbSLois Curfman McInnes 3065e76c431SBarry Smith Input Parameters: 307*c7afd0dbSLois Curfman McInnes + snes - nonlinear context 3085e76c431SBarry Smith . x - current iterate 3095e76c431SBarry Smith . f - residual evaluated at x 3105e76c431SBarry Smith . y - search direction (contains new iterate on output) 3115e76c431SBarry Smith . w - work vector 312*c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3135e76c431SBarry Smith 314393d2d9aSLois Curfman McInnes Output Parameters: 315*c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3165e76c431SBarry Smith . y - new iterate (contains search direction on input) 3175e76c431SBarry Smith . gnorm - 2-norm of g 3185e76c431SBarry Smith . ynorm - 2-norm of search length 319*c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 320fee21e36SBarry Smith 321c4a48953SLois Curfman McInnes Options Database Key: 322*c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 323c4a48953SLois Curfman McInnes 3245e76c431SBarry Smith Notes: 3255e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3265e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3275e76c431SBarry Smith 32828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 32928ae5a21SLois Curfman McInnes 33077c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 3315e76c431SBarry Smith @*/ 3325e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w, 333761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3345e76c431SBarry Smith { 335406556e6SLois Curfman McInnes /* 336406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 337406556e6SLois Curfman McInnes minimization problem: 338406556e6SLois Curfman McInnes min z(x): R^n -> R, 339406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 340406556e6SLois Curfman McInnes */ 341406556e6SLois Curfman McInnes 342ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 343ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 3443a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 3455e42470aSBarry Smith Scalar cinitslope, clambda; 3466b5873e3SBarry Smith #endif 3475e42470aSBarry Smith int ierr, count; 3485e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3495334005bSBarry Smith Scalar mone = -1.0,scale; 3505e76c431SBarry Smith 3513a40ed3dSBarry Smith PetscFunctionBegin; 3527857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 353761aaf1bSLois Curfman McInnes *flag = 0; 3545e76c431SBarry Smith alpha = neP->alpha; 3555e76c431SBarry Smith maxstep = neP->maxstep; 3565e76c431SBarry Smith steptol = neP->steptol; 3575e76c431SBarry Smith 358cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); 359a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 360a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 361a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 362a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 363a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 364ad922ac9SBarry Smith goto theend1; 36594a9d846SBarry Smith } 3665e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3675e42470aSBarry Smith scale = maxstep/(*ynorm); 3683a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 36994a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale)); 3706b5873e3SBarry Smith #else 37194a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 3726b5873e3SBarry Smith #endif 373393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 3745e76c431SBarry Smith *ynorm = maxstep; 3755e76c431SBarry Smith } 3765e76c431SBarry Smith minlambda = steptol/(*ynorm); 377a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 3783a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 379a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 3805e42470aSBarry Smith initslope = real(cinitslope); 3815e42470aSBarry Smith #else 382a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 3835e42470aSBarry Smith #endif 3845e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 3855e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 3865e76c431SBarry Smith 387393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 3885334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 38978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 390cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); 391406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 392393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 39394a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 394ad922ac9SBarry Smith goto theend1; 3955e76c431SBarry Smith } 3965e76c431SBarry Smith 3975e76c431SBarry Smith /* Fit points with quadratic */ 3985e76c431SBarry Smith lambda = 1.0; count = 0; 399a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4005e76c431SBarry Smith lambdaprev = lambda; 4015e76c431SBarry Smith gnormprev = *gnorm; 402ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 403ddd12767SBarry Smith else lambda = lambdatemp; 404393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 405ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4063a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 407ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4085e42470aSBarry Smith #else 409ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4105e42470aSBarry Smith #endif 41178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 412cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 413406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 414393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 415ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 416ad922ac9SBarry Smith goto theend1; 4175e76c431SBarry Smith } 4185e76c431SBarry Smith 4195e76c431SBarry Smith /* Fit points with cubic */ 4205e76c431SBarry Smith count = 1; 4215e76c431SBarry Smith while (1) { 4225e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4232b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4242b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 425ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 426393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 427761aaf1bSLois Curfman McInnes *flag = -1; break; 4285e76c431SBarry Smith } 429406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 430406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 431ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4322b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4335e76c431SBarry Smith d = b*b - 3*a*initslope; 4345e76c431SBarry Smith if (d < 0.0) d = 0.0; 4355e76c431SBarry Smith if (a == 0.0) { 4365e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4375e76c431SBarry Smith } else { 4385e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4395e76c431SBarry Smith } 4405e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4415e76c431SBarry Smith lambdatemp = .5*lambda; 4425e76c431SBarry Smith } 4435e76c431SBarry Smith lambdaprev = lambda; 4445e76c431SBarry Smith gnormprev = *gnorm; 4455e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4465e76c431SBarry Smith lambda = .1*lambda; 4475e76c431SBarry Smith } 4485e76c431SBarry Smith else lambda = lambdatemp; 449393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 450ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4513a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 452ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 453393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4545e42470aSBarry Smith #else 455ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4565e42470aSBarry Smith #endif 45778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 458cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 459406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 460393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 461ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4622715565aSLois Curfman McInnes goto theend1; 4632b022350SLois Curfman McInnes } else { 4642b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 4655e76c431SBarry Smith } 4665e76c431SBarry Smith count++; 4675e76c431SBarry Smith } 468ad922ac9SBarry Smith theend1: 4697857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 4715e76c431SBarry Smith } 47204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4735615d1e5SSatish Balay #undef __FUNC__ 4745615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 4754b828684SBarry Smith /*@C 476f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 4775e76c431SBarry Smith 478*c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 479*c7afd0dbSLois Curfman McInnes 4805e42470aSBarry Smith Input Parameters: 481*c7afd0dbSLois Curfman McInnes + snes - the SNES context 4825e42470aSBarry Smith . x - current iterate 4835e42470aSBarry Smith . f - residual evaluated at x 4845e42470aSBarry Smith . y - search direction (contains new iterate on output) 4855e42470aSBarry Smith . w - work vector 486*c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4875e42470aSBarry Smith 488c4a48953SLois Curfman McInnes Output Parameters: 489*c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4905e42470aSBarry Smith . y - new iterate (contains search direction on input) 4915e42470aSBarry Smith . gnorm - 2-norm of g 4925e42470aSBarry Smith . ynorm - 2-norm of search length 493*c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 494fee21e36SBarry Smith 495c4a48953SLois Curfman McInnes Options Database Key: 496*c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 497c4a48953SLois Curfman McInnes 4985e42470aSBarry Smith Notes: 499*c7afd0dbSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method. 5005e42470aSBarry Smith 501f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 502f59ffdedSLois Curfman McInnes 50377c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 5045e42470aSBarry Smith @*/ 5055e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 506761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5075e76c431SBarry Smith { 508406556e6SLois Curfman McInnes /* 509406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 510406556e6SLois Curfman McInnes minimization problem: 511406556e6SLois Curfman McInnes min z(x): R^n -> R, 512406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 513406556e6SLois Curfman McInnes */ 51440011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 5153a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 5165e42470aSBarry Smith Scalar cinitslope,clambda; 5176b5873e3SBarry Smith #endif 5185e42470aSBarry Smith int ierr,count; 5195e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 5205334005bSBarry Smith Scalar mone = -1.0,scale; 5215e76c431SBarry Smith 5223a40ed3dSBarry Smith PetscFunctionBegin; 5237857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 524761aaf1bSLois Curfman McInnes *flag = 0; 5255e42470aSBarry Smith alpha = neP->alpha; 5265e42470aSBarry Smith maxstep = neP->maxstep; 5275e42470aSBarry Smith steptol = neP->steptol; 5285e76c431SBarry Smith 529cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 530b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 53194a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 532b37302c6SLois Curfman McInnes *gnorm = fnorm; 533b37302c6SLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 534b37302c6SLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 535ad922ac9SBarry Smith goto theend2; 53694a9d846SBarry Smith } 5375e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5385e42470aSBarry Smith scale = maxstep/(*ynorm); 539393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 5405e42470aSBarry Smith *ynorm = maxstep; 5415e76c431SBarry Smith } 5425e42470aSBarry Smith minlambda = steptol/(*ynorm); 543a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 5443a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 545a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 5465e42470aSBarry Smith initslope = real(cinitslope); 5475e42470aSBarry Smith #else 548a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 5495e42470aSBarry Smith #endif 5505e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 5515e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 5525e42470aSBarry Smith 553393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 5545334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 55578b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 556cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 557406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 558393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 55994a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 560ad922ac9SBarry Smith goto theend2; 5615e42470aSBarry Smith } 5625e42470aSBarry Smith 5635e42470aSBarry Smith /* Fit points with quadratic */ 5645e42470aSBarry Smith lambda = 1.0; count = 0; 5655e42470aSBarry Smith count = 1; 5665e42470aSBarry Smith while (1) { 5675e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 568981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 569981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 570ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 571393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 572761aaf1bSLois Curfman McInnes *flag = -1; break; 5735e42470aSBarry Smith } 574a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5755e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 5765e42470aSBarry Smith lambda = .1*lambda; 5775e42470aSBarry Smith } else lambda = lambdatemp; 578393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 5795334005bSBarry Smith lambda = -lambda; 5803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 581393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 5825e42470aSBarry Smith #else 583393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 5845e42470aSBarry Smith #endif 58578b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 586cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 587406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 588393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 589981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 59006259719SBarry Smith break; 5915e42470aSBarry Smith } 5925e42470aSBarry Smith count++; 5935e42470aSBarry Smith } 594ad922ac9SBarry Smith theend2: 5957857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5963a40ed3dSBarry Smith PetscFunctionReturn(0); 5975e76c431SBarry Smith } 59804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5995615d1e5SSatish Balay #undef __FUNC__ 600d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 601c9e6a524SLois Curfman McInnes /*@C 60277c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 603f63b844aSLois Curfman McInnes by the method SNES_EQ_LS. 6045e76c431SBarry Smith 6055e76c431SBarry Smith Input Parameters: 606*c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 607*c7afd0dbSLois Curfman McInnes - func - pointer to int function 6085e76c431SBarry Smith 609fee21e36SBarry Smith Collective on SNES 610fee21e36SBarry Smith 611c4a48953SLois Curfman McInnes Available Routines: 612*c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 613f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 614f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 615*c7afd0dbSLois Curfman McInnes - SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel) 6165e76c431SBarry Smith 617c4a48953SLois Curfman McInnes Options Database Keys: 618*c7afd0dbSLois Curfman McInnes + -snes_eq_ls [basic,quadratic,cubic] - Selects line search 619*c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 620*c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 621*c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 622c4a48953SLois Curfman McInnes 6235e76c431SBarry Smith Calling sequence of func: 624f59ffdedSLois Curfman McInnes func (SNES snes, Vec x, Vec f, Vec g, Vec y, 625761aaf1bSLois Curfman McInnes Vec w, double fnorm, double *ynorm, 626761aaf1bSLois Curfman McInnes double *gnorm, *flag) 6275e76c431SBarry Smith 6285e76c431SBarry Smith Input parameters for func: 629*c7afd0dbSLois Curfman McInnes + snes - nonlinear context 6305e76c431SBarry Smith . x - current iterate 6315e76c431SBarry Smith . f - residual evaluated at x 6325e76c431SBarry Smith . y - search direction (contains new iterate on output) 6335e76c431SBarry Smith . w - work vector 634*c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6355e76c431SBarry Smith 6365e76c431SBarry Smith Output parameters for func: 637*c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6385e76c431SBarry Smith . y - new iterate (contains search direction on input) 6395e76c431SBarry Smith . gnorm - 2-norm of g 6405e76c431SBarry Smith . ynorm - 2-norm of search length 641*c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 642761aaf1bSLois Curfman McInnes on failure. 643f59ffdedSLois Curfman McInnes 644f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 645f59ffdedSLois Curfman McInnes 646f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch() 6475e76c431SBarry Smith @*/ 64877c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 649761aaf1bSLois Curfman McInnes double,double*,double*,int*)) 6505e76c431SBarry Smith { 65182bf6240SBarry Smith int ierr, (*f)(SNES,int (*f)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*)); 65282bf6240SBarry Smith 6533a40ed3dSBarry Smith PetscFunctionBegin; 654e1311b90SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch",(void **)&f);CHKERRQ(ierr); 65582bf6240SBarry Smith if (f) { 65682bf6240SBarry Smith ierr = (*f)(snes,func);CHKERRQ(ierr); 65782bf6240SBarry Smith } 6583a40ed3dSBarry Smith PetscFunctionReturn(0); 6595e76c431SBarry Smith } 66004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 66182bf6240SBarry Smith #undef __FUNC__ 66282bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS" 66382bf6240SBarry Smith int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 66482bf6240SBarry Smith double,double*,double*,int*)) 66582bf6240SBarry Smith { 66682bf6240SBarry Smith PetscFunctionBegin; 66782bf6240SBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 66882bf6240SBarry Smith PetscFunctionReturn(0); 66982bf6240SBarry Smith } 67004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 67104d965bbSLois Curfman McInnes /* 67204d965bbSLois Curfman McInnes SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 67382bf6240SBarry Smith 67404d965bbSLois Curfman McInnes Input Parameter: 67504d965bbSLois Curfman McInnes . snes - the SNES context 67604d965bbSLois Curfman McInnes 67704d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 67804d965bbSLois Curfman McInnes */ 67904d965bbSLois Curfman McInnes 6805615d1e5SSatish Balay #undef __FUNC__ 681d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS" 682f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 6835e42470aSBarry Smith { 6845e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 6856daaf66cSBarry Smith 6863a40ed3dSBarry Smith PetscFunctionBegin; 68776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n"); 68876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [basic,quadratic,cubic]\n",p); 68976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 69076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 69176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol); 6923a40ed3dSBarry Smith PetscFunctionReturn(0); 6935e42470aSBarry Smith } 69404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 69504d965bbSLois Curfman McInnes /* 69604d965bbSLois Curfman McInnes SNESView_EQ_LS - Prints the SNES_EQ_LS data structure. 69704d965bbSLois Curfman McInnes 69804d965bbSLois Curfman McInnes Input Parameters: 69904d965bbSLois Curfman McInnes . SNES - the SNES context 70004d965bbSLois Curfman McInnes . viewer - visualization context 70104d965bbSLois Curfman McInnes 70204d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 70304d965bbSLois Curfman McInnes */ 7045615d1e5SSatish Balay #undef __FUNC__ 705d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS" 706e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 707a935fc98SLois Curfman McInnes { 708a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 709d636dbe3SBarry Smith FILE *fd; 71019bcc07fSBarry Smith char *cstr; 71151695f54SLois Curfman McInnes int ierr; 71219bcc07fSBarry Smith ViewerType vtype; 713a935fc98SLois Curfman McInnes 7143a40ed3dSBarry Smith PetscFunctionBegin; 71519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 71619bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 71790ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 71819bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 71919bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 72019bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 72119bcc07fSBarry Smith else cstr = "unknown"; 72277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," line search variant: %s\n",cstr); 72377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," alpha=%g, maxstep=%g, steptol=%g\n", 724a935fc98SLois Curfman McInnes ls->alpha,ls->maxstep,ls->steptol); 7255cd90555SBarry Smith } else { 7265cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 72719bcc07fSBarry Smith } 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 729a935fc98SLois Curfman McInnes } 73004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 73104d965bbSLois Curfman McInnes /* 73204d965bbSLois Curfman McInnes SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method. 73304d965bbSLois Curfman McInnes 73404d965bbSLois Curfman McInnes Input Parameter: 73504d965bbSLois Curfman McInnes . snes - the SNES context 73604d965bbSLois Curfman McInnes 73704d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 73804d965bbSLois Curfman McInnes */ 7395615d1e5SSatish Balay #undef __FUNC__ 7405615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS" 741f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 7425e42470aSBarry Smith { 7435e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 7445e42470aSBarry Smith char ver[16]; 7455e42470aSBarry Smith double tmp; 74625cdf11fSBarry Smith int ierr,flg; 7475e42470aSBarry Smith 7483a40ed3dSBarry Smith PetscFunctionBegin; 74909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr); 75025cdf11fSBarry Smith if (flg) { 7515e42470aSBarry Smith ls->alpha = tmp; 7525e42470aSBarry Smith } 75309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr); 75425cdf11fSBarry Smith if (flg) { 7555e42470aSBarry Smith ls->maxstep = tmp; 7565e42470aSBarry Smith } 75709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr); 75825cdf11fSBarry Smith if (flg) { 7595e42470aSBarry Smith ls->steptol = tmp; 7605e42470aSBarry Smith } 76109d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr); 76225cdf11fSBarry Smith if (flg) { 76348d91487SBarry Smith if (!PetscStrcmp(ver,"basic")) { 76477c4ece6SBarry Smith SNESSetLineSearch(snes,SNESNoLineSearch); 7655e42470aSBarry Smith } 76629e0b56fSBarry Smith else if (!PetscStrcmp(ver,"basicnonorms")) { 76729e0b56fSBarry Smith SNESSetLineSearch(snes,SNESNoLineSearchNoNorms); 76829e0b56fSBarry Smith } 76948d91487SBarry Smith else if (!PetscStrcmp(ver,"quadratic")) { 77077c4ece6SBarry Smith SNESSetLineSearch(snes,SNESQuadraticLineSearch); 7715e42470aSBarry Smith } 77248d91487SBarry Smith else if (!PetscStrcmp(ver,"cubic")) { 77377c4ece6SBarry Smith SNESSetLineSearch(snes,SNESCubicLineSearch); 7745e42470aSBarry Smith } 775a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 7765e42470aSBarry Smith } 7773a40ed3dSBarry Smith PetscFunctionReturn(0); 7785e42470aSBarry Smith } 77904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 78004d965bbSLois Curfman McInnes /* 78104d965bbSLois Curfman McInnes SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method, 78204d965bbSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 78304d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 78404d965bbSLois Curfman McInnes 78504d965bbSLois Curfman McInnes 78604d965bbSLois Curfman McInnes Input Parameter: 78704d965bbSLois Curfman McInnes . snes - the SNES context 78804d965bbSLois Curfman McInnes 78904d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 79004d965bbSLois Curfman McInnes */ 7915615d1e5SSatish Balay #undef __FUNC__ 7925615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 793f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 7945e42470aSBarry Smith { 79582bf6240SBarry Smith int ierr; 7965e42470aSBarry Smith SNES_LS *neP; 7975e42470aSBarry Smith 7983a40ed3dSBarry Smith PetscFunctionBegin; 799a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 800a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 801a8c6a408SBarry Smith } 80282bf6240SBarry Smith 803f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 804f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 805f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 806f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 807f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 808f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 809f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 8105baf8537SBarry Smith snes->nwork = 0; 8115e42470aSBarry Smith 8120452661fSBarry Smith neP = PetscNew(SNES_LS); CHKPTRQ(neP); 813ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 8145e42470aSBarry Smith snes->data = (void *) neP; 8155e42470aSBarry Smith neP->alpha = 1.e-4; 8165e42470aSBarry Smith neP->maxstep = 1.e8; 8175e42470aSBarry Smith neP->steptol = 1.e-12; 8185e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 81982bf6240SBarry Smith 820e1311b90SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch","SNESSetLineSearch_LS", 821e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 82282bf6240SBarry Smith 8233a40ed3dSBarry Smith PetscFunctionReturn(0); 8245e42470aSBarry Smith } 8255e42470aSBarry Smith 82682bf6240SBarry Smith 82782bf6240SBarry Smith 82882bf6240SBarry Smith 829