1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3f1db9ecSBarry Smith static char vcid[] = "$Id: ls.c,v 1.117 1998/12/03 04:05:32 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; 69761aaf1bSLois Curfman McInnes int maxits, i, history_len, ierr, lits, lsfail; 70112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 716b5873e3SBarry Smith double fnorm, gnorm, xnorm, ynorm, *history; 725e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 735e76c431SBarry Smith 743a40ed3dSBarry Smith PetscFunctionBegin; 755e42470aSBarry Smith history = snes->conv_hist; /* convergence history */ 7651979daaSLois Curfman McInnes history_len = snes->conv_hist_size; /* convergence history length */ 775e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 785e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 7939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 805e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 815e42470aSBarry Smith G = snes->work[1]; 825e42470aSBarry Smith W = snes->work[2]; 835e76c431SBarry Smith 8425ed9cd1SBarry Smith snes->iter = 0; 855334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr); /* F(X) */ 86cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr); /* fnorm <- ||F|| */ 875e42470aSBarry Smith snes->norm = fnorm; 8851979daaSLois Curfman McInnes if (history) history[0] = fnorm; 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++ ) { 975e42470aSBarry Smith snes->iter = i+1; 985e76c431SBarry Smith 99ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1005334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr); 1015334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr); 10278b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr); 1037a00f4a9SLois Curfman McInnes snes->linear_its += PetscAbsInt(lits); 104af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 105ea4d3ed3SLois Curfman McInnes 106ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 107ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 108ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 109ea4d3ed3SLois Curfman McInnes */ 11081b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr); 111ddd12767SBarry Smith ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr); 112af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 113761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 1145e76c431SBarry Smith 11539e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 11639e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1175e76c431SBarry Smith fnorm = gnorm; 1185e76c431SBarry Smith 1195e42470aSBarry Smith snes->norm = fnorm; 1205e76c431SBarry Smith if (history && history_len > i+1) history[i+1] = fnorm; 12194a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1225e76c431SBarry Smith 1235e76c431SBarry Smith /* Test for convergence */ 12429e0b56fSBarry Smith if (snes->converged) { 12529e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr); /* xnorm = || X || */ 126bbb6d6a8SBarry Smith if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) { 12716c95cacSBarry Smith break; 12816c95cacSBarry Smith } 12916c95cacSBarry Smith } 13029e0b56fSBarry Smith } 13139e2f89bSBarry Smith if (X != snes->vec_sol) { 132393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr); 13339e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 13439e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 13539e2f89bSBarry Smith } 13652392280SLois Curfman McInnes if (i == maxits) { 137981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 13852392280SLois Curfman McInnes i--; 13952392280SLois Curfman McInnes } 14051979daaSLois Curfman McInnes if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1; 1415e42470aSBarry Smith *outits = i+1; 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 1435e76c431SBarry Smith } 14404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 14504d965bbSLois Curfman McInnes /* 14604d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 14704d965bbSLois Curfman McInnes of the SNES_EQ_LS nonlinear solver. 14804d965bbSLois Curfman McInnes 14904d965bbSLois Curfman McInnes Input Parameter: 15004d965bbSLois Curfman McInnes . snes - the SNES context 15104d965bbSLois Curfman McInnes . x - the solution vector 15204d965bbSLois Curfman McInnes 15304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 15404d965bbSLois Curfman McInnes 15504d965bbSLois Curfman McInnes Notes: 15604d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 15704d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 15804d965bbSLois Curfman McInnes the call to SNESSolve(). 15904d965bbSLois Curfman McInnes */ 1605615d1e5SSatish Balay #undef __FUNC__ 1615615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS" 162f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1635e76c431SBarry Smith { 1645e42470aSBarry Smith int ierr; 1653a40ed3dSBarry Smith 1663a40ed3dSBarry Smith PetscFunctionBegin; 16781b6cf68SLois Curfman McInnes snes->nwork = 4; 168d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1695e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 17081b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1713a40ed3dSBarry Smith PetscFunctionReturn(0); 1725e76c431SBarry Smith } 17304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 17404d965bbSLois Curfman McInnes /* 17504d965bbSLois Curfman McInnes SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created 17604d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 17704d965bbSLois Curfman McInnes 17804d965bbSLois Curfman McInnes Input Parameter: 17904d965bbSLois Curfman McInnes . snes - the SNES context 18004d965bbSLois Curfman McInnes 181de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 18204d965bbSLois Curfman McInnes */ 1835615d1e5SSatish Balay #undef __FUNC__ 184d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS" 185e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 1865e76c431SBarry Smith { 187393d2d9aSLois Curfman McInnes int ierr; 1883a40ed3dSBarry Smith 1893a40ed3dSBarry Smith PetscFunctionBegin; 1905baf8537SBarry Smith if (snes->nwork) { 1914b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr); 19221c89e3eSBarry Smith } 1930452661fSBarry Smith PetscFree(snes->data); 1943a40ed3dSBarry Smith PetscFunctionReturn(0); 1955e76c431SBarry Smith } 19604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1975615d1e5SSatish Balay #undef __FUNC__ 1985615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch" 19982bf6240SBarry Smith 2004b828684SBarry Smith /*@C 2015e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2025e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2035e76c431SBarry Smith to serve as a template and is not recommended for general use. 2045e76c431SBarry Smith 205c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 206c7afd0dbSLois Curfman McInnes 2075e76c431SBarry Smith Input Parameters: 208c7afd0dbSLois Curfman McInnes + snes - nonlinear context 2095e76c431SBarry Smith . x - current iterate 2105e76c431SBarry Smith . f - residual evaluated at x 2115e76c431SBarry Smith . y - search direction (contains new iterate on output) 2125e76c431SBarry Smith . w - work vector 213c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2145e76c431SBarry Smith 215c4a48953SLois Curfman McInnes Output Parameters: 216c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2175e76c431SBarry Smith . y - new iterate (contains search direction on input) 2185e76c431SBarry Smith . gnorm - 2-norm of g 2195e76c431SBarry Smith . ynorm - 2-norm of search length 220c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 221fee21e36SBarry Smith 222c4a48953SLois Curfman McInnes Options Database Key: 223c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 224c4a48953SLois Curfman McInnes 22528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 22628ae5a21SLois Curfman McInnes 227f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 22877c4ece6SBarry Smith SNESSetLineSearch() 2295e76c431SBarry Smith @*/ 2305e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 231761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 2325e76c431SBarry Smith { 2335e42470aSBarry Smith int ierr; 2345334005bSBarry Smith Scalar mone = -1.0; 2355334005bSBarry Smith 2363a40ed3dSBarry Smith PetscFunctionBegin; 237761aaf1bSLois Curfman McInnes *flag = 0; 2387857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 239cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 240ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 241ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 242cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 2437857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2455e76c431SBarry Smith } 24604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2475615d1e5SSatish Balay #undef __FUNC__ 24829e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms" 24982bf6240SBarry Smith 25029e0b56fSBarry Smith /*@C 25129e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 25229e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 25329e0b56fSBarry Smith even compute the norm of the function or search direction; this 25429e0b56fSBarry Smith is intended only when you know the full step is fine and are 25529e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 256c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 257c7afd0dbSLois Curfman McInnes 258c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 25929e0b56fSBarry Smith 26029e0b56fSBarry Smith Input Parameters: 261c7afd0dbSLois Curfman McInnes + snes - nonlinear context 26229e0b56fSBarry Smith . x - current iterate 26329e0b56fSBarry Smith . f - residual evaluated at x 26429e0b56fSBarry Smith . y - search direction (contains new iterate on output) 26529e0b56fSBarry Smith . w - work vector 266c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 26729e0b56fSBarry Smith 26829e0b56fSBarry Smith Output Parameters: 269c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 27029e0b56fSBarry Smith . gnorm - not changed 27129e0b56fSBarry Smith . ynorm - not changed 272c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 273fee21e36SBarry Smith 27429e0b56fSBarry Smith Options Database Key: 275c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 27629e0b56fSBarry Smith 27729e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 27829e0b56fSBarry Smith 27929e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 28029e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 28129e0b56fSBarry Smith @*/ 28229e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 28329e0b56fSBarry Smith double fnorm, double *ynorm, double *gnorm,int *flag ) 28429e0b56fSBarry Smith { 28529e0b56fSBarry Smith int ierr; 28629e0b56fSBarry Smith Scalar mone = -1.0; 28729e0b56fSBarry Smith 2883a40ed3dSBarry Smith PetscFunctionBegin; 28929e0b56fSBarry Smith *flag = 0; 29029e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 29129e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr); /* y <- y - x */ 29229e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y) */ 29329e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 29529e0b56fSBarry Smith } 29604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 29729e0b56fSBarry Smith #undef __FUNC__ 2985615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch" 2994b828684SBarry Smith /*@C 300f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3015e76c431SBarry Smith 302c7afd0dbSLois Curfman McInnes Collective on SNES 303c7afd0dbSLois Curfman McInnes 3045e76c431SBarry Smith Input Parameters: 305c7afd0dbSLois Curfman McInnes + snes - nonlinear context 3065e76c431SBarry Smith . x - current iterate 3075e76c431SBarry Smith . f - residual evaluated at x 3085e76c431SBarry Smith . y - search direction (contains new iterate on output) 3095e76c431SBarry Smith . w - work vector 310c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3115e76c431SBarry Smith 312393d2d9aSLois Curfman McInnes Output Parameters: 313c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3145e76c431SBarry Smith . y - new iterate (contains search direction on input) 3155e76c431SBarry Smith . gnorm - 2-norm of g 3165e76c431SBarry Smith . ynorm - 2-norm of search length 317c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 318fee21e36SBarry Smith 319c4a48953SLois Curfman McInnes Options Database Key: 320c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 321c4a48953SLois Curfman McInnes 3225e76c431SBarry Smith Notes: 3235e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3245e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3255e76c431SBarry Smith 32628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 32728ae5a21SLois Curfman McInnes 32877c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 3295e76c431SBarry Smith @*/ 3305e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w, 331761aaf1bSLois Curfman McInnes double fnorm,double *ynorm,double *gnorm,int *flag) 3325e76c431SBarry Smith { 333406556e6SLois Curfman McInnes /* 334406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 335406556e6SLois Curfman McInnes minimization problem: 336406556e6SLois Curfman McInnes min z(x): R^n -> R, 337406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 338406556e6SLois Curfman McInnes */ 339406556e6SLois Curfman McInnes 340ddd12767SBarry Smith double steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2; 341ea4d3ed3SLois Curfman McInnes double maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg; 3423a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 3435e42470aSBarry Smith Scalar cinitslope, clambda; 3446b5873e3SBarry Smith #endif 3455e42470aSBarry Smith int ierr, count; 3465e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3475334005bSBarry Smith Scalar mone = -1.0,scale; 3485e76c431SBarry Smith 3493a40ed3dSBarry Smith PetscFunctionBegin; 3507857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 351761aaf1bSLois Curfman McInnes *flag = 0; 3525e76c431SBarry Smith alpha = neP->alpha; 3535e76c431SBarry Smith maxstep = neP->maxstep; 3545e76c431SBarry Smith steptol = neP->steptol; 3555e76c431SBarry Smith 356cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr); 357a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 358a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 359a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 360a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 361a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 362ad922ac9SBarry Smith goto theend1; 36394a9d846SBarry Smith } 3645e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3655e42470aSBarry Smith scale = maxstep/(*ynorm); 3663a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 36792318cfeSSatish Balay PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale)); 3686b5873e3SBarry Smith #else 36994a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 3706b5873e3SBarry Smith #endif 371393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 3725e76c431SBarry Smith *ynorm = maxstep; 3735e76c431SBarry Smith } 3745e76c431SBarry Smith minlambda = steptol/(*ynorm); 375a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 3763a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 377a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 37892318cfeSSatish Balay initslope = PetscReal(cinitslope); 3795e42470aSBarry Smith #else 380a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 3815e42470aSBarry Smith #endif 3825e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 3835e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 3845e76c431SBarry Smith 385393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 3865334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 38778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 388cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); 389406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 390393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 39194a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 392ad922ac9SBarry Smith goto theend1; 3935e76c431SBarry Smith } 3945e76c431SBarry Smith 3955e76c431SBarry Smith /* Fit points with quadratic */ 3965e76c431SBarry Smith lambda = 1.0; count = 0; 397a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 3985e76c431SBarry Smith lambdaprev = lambda; 3995e76c431SBarry Smith gnormprev = *gnorm; 400ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 401ddd12767SBarry Smith else lambda = lambdatemp; 402393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 403ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4043a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 405ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4065e42470aSBarry Smith #else 407ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4085e42470aSBarry Smith #endif 40978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 410cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 411406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 412393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 413ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 414ad922ac9SBarry Smith goto theend1; 4155e76c431SBarry Smith } 4165e76c431SBarry Smith 4175e76c431SBarry Smith /* Fit points with cubic */ 4185e76c431SBarry Smith count = 1; 4195e76c431SBarry Smith while (1) { 4205e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4212b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 4222b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 423ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 424393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 425761aaf1bSLois Curfman McInnes *flag = -1; break; 4265e76c431SBarry Smith } 427406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 428406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 429ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4302b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4315e76c431SBarry Smith d = b*b - 3*a*initslope; 4325e76c431SBarry Smith if (d < 0.0) d = 0.0; 4335e76c431SBarry Smith if (a == 0.0) { 4345e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4355e76c431SBarry Smith } else { 4365e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4375e76c431SBarry Smith } 4385e76c431SBarry Smith if (lambdatemp > .5*lambda) { 4395e76c431SBarry Smith lambdatemp = .5*lambda; 4405e76c431SBarry Smith } 4415e76c431SBarry Smith lambdaprev = lambda; 4425e76c431SBarry Smith gnormprev = *gnorm; 4435e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 4445e76c431SBarry Smith lambda = .1*lambda; 4455e76c431SBarry Smith } 4465e76c431SBarry Smith else lambda = lambdatemp; 447393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 448ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 4493a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 450ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 451393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4525e42470aSBarry Smith #else 453ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr); 4545e42470aSBarry Smith #endif 45578b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 456cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 457406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 458393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 459ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 4602715565aSLois Curfman McInnes goto theend1; 4612b022350SLois Curfman McInnes } else { 4622b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 4635e76c431SBarry Smith } 4645e76c431SBarry Smith count++; 4655e76c431SBarry Smith } 466ad922ac9SBarry Smith theend1: 4677857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 4695e76c431SBarry Smith } 47004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4715615d1e5SSatish Balay #undef __FUNC__ 4725615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch" 4734b828684SBarry Smith /*@C 474f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 4755e76c431SBarry Smith 476c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 477c7afd0dbSLois Curfman McInnes 4785e42470aSBarry Smith Input Parameters: 479c7afd0dbSLois Curfman McInnes + snes - the SNES context 4805e42470aSBarry Smith . x - current iterate 4815e42470aSBarry Smith . f - residual evaluated at x 4825e42470aSBarry Smith . y - search direction (contains new iterate on output) 4835e42470aSBarry Smith . w - work vector 484c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4855e42470aSBarry Smith 486c4a48953SLois Curfman McInnes Output Parameters: 487c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4885e42470aSBarry Smith . y - new iterate (contains search direction on input) 4895e42470aSBarry Smith . gnorm - 2-norm of g 4905e42470aSBarry Smith . ynorm - 2-norm of search length 491c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 492fee21e36SBarry Smith 493c4a48953SLois Curfman McInnes Options Database Key: 494c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 495c4a48953SLois Curfman McInnes 4965e42470aSBarry Smith Notes: 497c7afd0dbSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method. 4985e42470aSBarry Smith 499f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 500f59ffdedSLois Curfman McInnes 50177c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch() 5025e42470aSBarry Smith @*/ 5035e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 504761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 5055e76c431SBarry Smith { 506406556e6SLois Curfman McInnes /* 507406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 508406556e6SLois Curfman McInnes minimization problem: 509406556e6SLois Curfman McInnes min z(x): R^n -> R, 510406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 511406556e6SLois Curfman McInnes */ 51240011551SBarry Smith double steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp; 5133a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 5145e42470aSBarry Smith Scalar cinitslope,clambda; 5156b5873e3SBarry Smith #endif 5165e42470aSBarry Smith int ierr,count; 5175e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 5185334005bSBarry Smith Scalar mone = -1.0,scale; 5195e76c431SBarry Smith 5203a40ed3dSBarry Smith PetscFunctionBegin; 5217857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 522761aaf1bSLois Curfman McInnes *flag = 0; 5235e42470aSBarry Smith alpha = neP->alpha; 5245e42470aSBarry Smith maxstep = neP->maxstep; 5255e42470aSBarry Smith steptol = neP->steptol; 5265e76c431SBarry Smith 527cddf8d76SBarry Smith VecNorm(y, NORM_2,ynorm ); 528b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 52994a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 530b37302c6SLois Curfman McInnes *gnorm = fnorm; 531b37302c6SLois Curfman McInnes ierr = VecCopy(x,y); CHKERRQ(ierr); 532b37302c6SLois Curfman McInnes ierr = VecCopy(f,g); CHKERRQ(ierr); 533ad922ac9SBarry Smith goto theend2; 53494a9d846SBarry Smith } 5355e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5365e42470aSBarry Smith scale = maxstep/(*ynorm); 537393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 5385e42470aSBarry Smith *ynorm = maxstep; 5395e76c431SBarry Smith } 5405e42470aSBarry Smith minlambda = steptol/(*ynorm); 541a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 5423a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 543a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 54492318cfeSSatish Balay initslope = PetscReal(cinitslope); 5455e42470aSBarry Smith #else 546a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 5475e42470aSBarry Smith #endif 5485e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 5495e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 5505e42470aSBarry Smith 551393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 5525334005bSBarry Smith ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr); 55378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 554cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 555406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 556393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 55794a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 558ad922ac9SBarry Smith goto theend2; 5595e42470aSBarry Smith } 5605e42470aSBarry Smith 5615e42470aSBarry Smith /* Fit points with quadratic */ 5625e42470aSBarry Smith lambda = 1.0; count = 0; 5635e42470aSBarry Smith count = 1; 5645e42470aSBarry Smith while (1) { 5655e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 566981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 567981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 568ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 569393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 570761aaf1bSLois Curfman McInnes *flag = -1; break; 5715e42470aSBarry Smith } 572a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5735e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 5745e42470aSBarry Smith lambda = .1*lambda; 5755e42470aSBarry Smith } else lambda = lambdatemp; 576393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 5775334005bSBarry Smith lambda = -lambda; 5783a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 579393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 5805e42470aSBarry Smith #else 581393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 5825e42470aSBarry Smith #endif 58378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 584cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr); 585406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 586393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 587981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 58806259719SBarry Smith break; 5895e42470aSBarry Smith } 5905e42470aSBarry Smith count++; 5915e42470aSBarry Smith } 592ad922ac9SBarry Smith theend2: 5937857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5943a40ed3dSBarry Smith PetscFunctionReturn(0); 5955e76c431SBarry Smith } 59604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5975615d1e5SSatish Balay #undef __FUNC__ 598d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch" 599c9e6a524SLois Curfman McInnes /*@C 60077c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 601f63b844aSLois Curfman McInnes by the method SNES_EQ_LS. 6025e76c431SBarry Smith 6035e76c431SBarry Smith Input Parameters: 604c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 605c7afd0dbSLois Curfman McInnes - func - pointer to int function 6065e76c431SBarry Smith 607fee21e36SBarry Smith Collective on SNES 608fee21e36SBarry Smith 609c4a48953SLois Curfman McInnes Available Routines: 610c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 611f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 612f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 613c7afd0dbSLois Curfman McInnes - SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel) 6145e76c431SBarry Smith 615c4a48953SLois Curfman McInnes Options Database Keys: 616c7afd0dbSLois Curfman McInnes + -snes_eq_ls [basic,quadratic,cubic] - Selects line search 617c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 618c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 619c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 620c4a48953SLois Curfman McInnes 6215e76c431SBarry Smith Calling sequence of func: 622bd208895SLois Curfman McInnes .vb 623bd208895SLois Curfman McInnes func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 624bd208895SLois Curfman McInnes double fnorm, double *ynorm, double *gnorm, *flag) 625bd208895SLois Curfman McInnes .ve 6265e76c431SBarry Smith 6275e76c431SBarry Smith Input parameters for func: 628c7afd0dbSLois Curfman McInnes + snes - nonlinear context 6295e76c431SBarry Smith . x - current iterate 6305e76c431SBarry Smith . f - residual evaluated at x 6315e76c431SBarry Smith . y - search direction (contains new iterate on output) 6325e76c431SBarry Smith . w - work vector 633c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6345e76c431SBarry Smith 6355e76c431SBarry Smith Output parameters for func: 636c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6375e76c431SBarry Smith . y - new iterate (contains search direction on input) 6385e76c431SBarry Smith . gnorm - 2-norm of g 6395e76c431SBarry Smith . ynorm - 2-norm of search length 640c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 641761aaf1bSLois Curfman McInnes on failure. 642f59ffdedSLois Curfman McInnes 643f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 644f59ffdedSLois Curfman McInnes 645f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch() 6465e76c431SBarry Smith @*/ 64777c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 648761aaf1bSLois Curfman McInnes double,double*,double*,int*)) 6495e76c431SBarry Smith { 650415aef20SBarry Smith int ierr, (*f)(SNES,int (*)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*)); 65182bf6240SBarry Smith 6523a40ed3dSBarry Smith PetscFunctionBegin; 65394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 65482bf6240SBarry Smith if (f) { 65582bf6240SBarry Smith ierr = (*f)(snes,func);CHKERRQ(ierr); 65682bf6240SBarry Smith } 6573a40ed3dSBarry Smith PetscFunctionReturn(0); 6585e76c431SBarry Smith } 65904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 660fb2e594dSBarry Smith EXTERN_C_BEGIN 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 } 670fb2e594dSBarry Smith EXTERN_C_END 67104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 67204d965bbSLois Curfman McInnes /* 67304d965bbSLois Curfman McInnes SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 67482bf6240SBarry Smith 67504d965bbSLois Curfman McInnes Input Parameter: 67604d965bbSLois Curfman McInnes . snes - the SNES context 67704d965bbSLois Curfman McInnes 67804d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 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 /* 696de49b36dSLois Curfman McInnes SNESView_EQ_LS - Prints info from 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); 716*3f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_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")) { 764b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESNoLineSearch);CHKERRQ(ierr); 765b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"basicnonorms")) { 766b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);CHKERRQ(ierr); 767b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"quadratic")) { 768b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch);CHKERRQ(ierr); 769b2d88d52SBarry Smith } else if (!PetscStrcmp(ver,"cubic")) { 770b2d88d52SBarry Smith ierr = SNESSetLineSearch(snes,SNESCubicLineSearch);CHKERRQ(ierr); 7715e42470aSBarry Smith } 772a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 7735e42470aSBarry Smith } 7743a40ed3dSBarry Smith PetscFunctionReturn(0); 7755e42470aSBarry Smith } 77604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 77704d965bbSLois Curfman McInnes /* 77804d965bbSLois Curfman McInnes SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method, 77904d965bbSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 78004d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 78104d965bbSLois Curfman McInnes 78204d965bbSLois Curfman McInnes 78304d965bbSLois Curfman McInnes Input Parameter: 78404d965bbSLois Curfman McInnes . snes - the SNES context 78504d965bbSLois Curfman McInnes 78604d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 78704d965bbSLois Curfman McInnes */ 788fb2e594dSBarry Smith EXTERN_C_BEGIN 7895615d1e5SSatish Balay #undef __FUNC__ 7905615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS" 791f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 7925e42470aSBarry Smith { 79382bf6240SBarry Smith int ierr; 7945e42470aSBarry Smith SNES_LS *neP; 7955e42470aSBarry Smith 7963a40ed3dSBarry Smith PetscFunctionBegin; 797a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 798a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 799a8c6a408SBarry Smith } 80082bf6240SBarry Smith 801f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 802f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 803f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 804f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 805f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 806f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 807f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 8085baf8537SBarry Smith snes->nwork = 0; 8095e42470aSBarry Smith 8100452661fSBarry Smith neP = PetscNew(SNES_LS); CHKPTRQ(neP); 811ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 8125e42470aSBarry Smith snes->data = (void *) neP; 8135e42470aSBarry Smith neP->alpha = 1.e-4; 8145e42470aSBarry Smith neP->maxstep = 1.e8; 8155e42470aSBarry Smith neP->steptol = 1.e-12; 8165e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 81782bf6240SBarry Smith 81894d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 819e1311b90SBarry Smith (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 82082bf6240SBarry Smith 8213a40ed3dSBarry Smith PetscFunctionReturn(0); 8225e42470aSBarry Smith } 823fb2e594dSBarry Smith EXTERN_C_END 82482bf6240SBarry Smith 82582bf6240SBarry Smith 82682bf6240SBarry Smith 827