1*90cbd584SBarry Smith /*$Id: ls.c,v 1.155 2000/05/04 14:03:59 balay Exp bsmith $*/ 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 604d965bbSLois Curfman McInnes 704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 804d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 1004d965bbSLois Curfman McInnes respectively. 1104d965bbSLois Curfman McInnes 12fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 1304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 1404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 15fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 1604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 1704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 1804d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 1904d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 2004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 2104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 2204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 2304d965bbSLois Curfman McInnes for all nonlinear solvers. 2404d965bbSLois Curfman McInnes 2504d965bbSLois Curfman McInnes Another key routine is: 2604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 2704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 2804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 2904d965bbSLois Curfman McInnes SNESSolve() if necessary. 3004d965bbSLois Curfman McInnes 3104d965bbSLois Curfman McInnes Additional basic routines are: 3204d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 3304d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 3404d965bbSLois Curfman McInnes have actually been used. 3504d965bbSLois Curfman McInnes These are called by application codes via the interface routines 3604d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 3704d965bbSLois Curfman McInnes 3804d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 3904d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 4004d965bbSLois Curfman McInnes above description applies to these categories also. 4104d965bbSLois Curfman McInnes 4204d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 435e42470aSBarry Smith /* 4404d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 4504d965bbSLois Curfman McInnes method with a line search. 465e76c431SBarry Smith 4704d965bbSLois Curfman McInnes Input Parameters: 4804d965bbSLois Curfman McInnes . snes - the SNES context 495e76c431SBarry Smith 5004d965bbSLois Curfman McInnes Output Parameter: 51c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 5204d965bbSLois Curfman McInnes 5304d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 545e76c431SBarry Smith 555e76c431SBarry Smith Notes: 565e76c431SBarry Smith This implements essentially a truncated Newton method with a 575e76c431SBarry Smith line search. By default a cubic backtracking line search 585e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 595e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 60393d2d9aSLois Curfman McInnes and Schnabel. 615e76c431SBarry Smith */ 625615d1e5SSatish Balay #undef __FUNC__ 63b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 64f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 655e76c431SBarry Smith { 666831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 6784c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 68112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 69329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 705e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 715e76c431SBarry Smith 723a40ed3dSBarry Smith PetscFunctionBegin; 73184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 74184914b5SBarry Smith 755e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 765e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 7739e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 785e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 795e42470aSBarry Smith G = snes->work[1]; 805e42470aSBarry Smith W = snes->work[2]; 815e76c431SBarry Smith 824c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 834c49b128SBarry Smith snes->iter = 0; 844c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 855334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 86cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 870f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 885e42470aSBarry Smith snes->norm = fnorm; 890f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 9084c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 9194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 925e76c431SBarry Smith 93184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; 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++) { 995e76c431SBarry Smith 100ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1015334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1025334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 10378b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 104*90cbd584SBarry Smith /* should check what happened to the linear solve? */ 1053505fcc1SBarry Smith snes->linear_its += 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); 113af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,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); 1153505fcc1SBarry Smith if (lsfail) { 116*90cbd584SBarry Smith double a1; 1173505fcc1SBarry Smith snes->nfailures++; 1183505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 119*90cbd584SBarry Smith /* Compute || J^T b||/|| b || */ 120*90cbd584SBarry Smith ierr = MatMultTranspose(snes->jacobian,F,W);CHKERRQ(ierr); 121*90cbd584SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 122*90cbd584SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: || J^T b|| %g || b || %g || J^T b||/|| b || %g\n",a1,fnorm,a1/fnorm); 123*90cbd584SBarry Smith /* || J^T b || == 0 implies the linear system is inconsistent and also that F() has a local minimum at the point */ 124*90cbd584SBarry Smith if (a1/fnorm < 1.e-4) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1253505fcc1SBarry Smith break; 1263505fcc1SBarry Smith } 1275e76c431SBarry Smith 12839e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 12939e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1305e76c431SBarry Smith fnorm = gnorm; 1315e76c431SBarry Smith 1320f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 13384c569c9SBarry Smith snes->iter = i+1; 1345e42470aSBarry Smith snes->norm = fnorm; 1350f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 13684c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 13794a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1385e76c431SBarry Smith 1395e76c431SBarry Smith /* Test for convergence */ 14029e0b56fSBarry Smith if (snes->converged) { 14129e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1423505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1433505fcc1SBarry Smith if (snes->reason) { 14416c95cacSBarry Smith break; 14516c95cacSBarry Smith } 14616c95cacSBarry Smith } 14729e0b56fSBarry Smith } 14839e2f89bSBarry Smith if (X != snes->vec_sol) { 149393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 150e15f7bb5SBarry Smith } 151e15f7bb5SBarry Smith if (F != snes->vec_func) { 152e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 153e15f7bb5SBarry Smith } 15439e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 15539e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 15652392280SLois Curfman McInnes if (i == maxits) { 157981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 15852392280SLois Curfman McInnes i--; 1593505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 16052392280SLois Curfman McInnes } 1610f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1620f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1635e42470aSBarry Smith *outits = i+1; 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 1655e76c431SBarry Smith } 16604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 16704d965bbSLois Curfman McInnes /* 16804d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 1696831982aSBarry Smith of the SNESEQLS nonlinear solver. 17004d965bbSLois Curfman McInnes 17104d965bbSLois Curfman McInnes Input Parameter: 17204d965bbSLois Curfman McInnes . snes - the SNES context 17304d965bbSLois Curfman McInnes . x - the solution vector 17404d965bbSLois Curfman McInnes 17504d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 17604d965bbSLois Curfman McInnes 17704d965bbSLois Curfman McInnes Notes: 17804d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 17904d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 18004d965bbSLois Curfman McInnes the call to SNESSolve(). 18104d965bbSLois Curfman McInnes */ 1825615d1e5SSatish Balay #undef __FUNC__ 183b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 184f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 1855e76c431SBarry Smith { 1865e42470aSBarry Smith int ierr; 1873a40ed3dSBarry Smith 1883a40ed3dSBarry Smith PetscFunctionBegin; 18981b6cf68SLois Curfman McInnes snes->nwork = 4; 190d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 1915e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 19281b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 1945e76c431SBarry Smith } 19504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 19604d965bbSLois Curfman McInnes /* 1976831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 19804d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 19904d965bbSLois Curfman McInnes 20004d965bbSLois Curfman McInnes Input Parameter: 20104d965bbSLois Curfman McInnes . snes - the SNES context 20204d965bbSLois Curfman McInnes 203de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 20404d965bbSLois Curfman McInnes */ 2055615d1e5SSatish Balay #undef __FUNC__ 206b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 207e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2085e76c431SBarry Smith { 209393d2d9aSLois Curfman McInnes int ierr; 2103a40ed3dSBarry Smith 2113a40ed3dSBarry Smith PetscFunctionBegin; 2125baf8537SBarry Smith if (snes->nwork) { 2134b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 21421c89e3eSBarry Smith } 2155c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2163a40ed3dSBarry Smith PetscFunctionReturn(0); 2175e76c431SBarry Smith } 21804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2195615d1e5SSatish Balay #undef __FUNC__ 220b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 22182bf6240SBarry Smith 2224b828684SBarry Smith /*@C 2235e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2245e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2255e76c431SBarry Smith to serve as a template and is not recommended for general use. 2265e76c431SBarry Smith 227c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 228c7afd0dbSLois Curfman McInnes 2295e76c431SBarry Smith Input Parameters: 230c7afd0dbSLois Curfman McInnes + snes - nonlinear context 231af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2325e76c431SBarry Smith . x - current iterate 2335e76c431SBarry Smith . f - residual evaluated at x 2345e76c431SBarry Smith . y - search direction (contains new iterate on output) 2355e76c431SBarry Smith . w - work vector 236c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2375e76c431SBarry Smith 238c4a48953SLois Curfman McInnes Output Parameters: 239c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2405e76c431SBarry Smith . y - new iterate (contains search direction on input) 2415e76c431SBarry Smith . gnorm - 2-norm of g 2425e76c431SBarry Smith . ynorm - 2-norm of search length 243c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 244fee21e36SBarry Smith 245c4a48953SLois Curfman McInnes Options Database Key: 246c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 247c4a48953SLois Curfman McInnes 24836851e7fSLois Curfman McInnes Level: advanced 24936851e7fSLois Curfman McInnes 25028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 25128ae5a21SLois Curfman McInnes 252f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 253af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2545e76c431SBarry Smith @*/ 255af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 256329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 2575e76c431SBarry Smith { 2585e42470aSBarry Smith int ierr; 2595334005bSBarry Smith Scalar mone = -1.0; 2606831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 2618f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2625334005bSBarry Smith 2633a40ed3dSBarry Smith PetscFunctionBegin; 264761aaf1bSLois Curfman McInnes *flag = 0; 2657857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 266a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 267ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 2688f99978dSLois Curfman McInnes if (neP->CheckStep) { 2698f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 2708f99978dSLois Curfman McInnes } 271ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 272a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 2737857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 2743a40ed3dSBarry Smith PetscFunctionReturn(0); 2755e76c431SBarry Smith } 27604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2775615d1e5SSatish Balay #undef __FUNC__ 278b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 27982bf6240SBarry Smith 28029e0b56fSBarry Smith /*@C 28129e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 28229e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 28329e0b56fSBarry Smith even compute the norm of the function or search direction; this 28429e0b56fSBarry Smith is intended only when you know the full step is fine and are 28529e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 286c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 287c7afd0dbSLois Curfman McInnes 288c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 28929e0b56fSBarry Smith 29029e0b56fSBarry Smith Input Parameters: 291c7afd0dbSLois Curfman McInnes + snes - nonlinear context 292af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 29329e0b56fSBarry Smith . x - current iterate 29429e0b56fSBarry Smith . f - residual evaluated at x 29529e0b56fSBarry Smith . y - search direction (contains new iterate on output) 29629e0b56fSBarry Smith . w - work vector 297c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 29829e0b56fSBarry Smith 29929e0b56fSBarry Smith Output Parameters: 300c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 30129e0b56fSBarry Smith . gnorm - not changed 30229e0b56fSBarry Smith . ynorm - not changed 303c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 304fee21e36SBarry Smith 30529e0b56fSBarry Smith Options Database Key: 306c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 30729e0b56fSBarry Smith 3088cbba510SBarry Smith Notes: 309ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 310ea56f5baSLois Curfman McInnes either the options 311ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 312ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 313329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 314329f5518SBarry Smith otherwise, the SNES solver will generate an error. 315329f5518SBarry Smith 316329f5518SBarry Smith During the final iteration this will not evaluate the function at 317329f5518SBarry Smith the solution point. This is to save a function evaluation while 318329f5518SBarry Smith using pseudo-timestepping. 3198cbba510SBarry Smith 320ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 321ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 322ea56f5baSLois Curfman McInnes correct, since they are not computed. 323ea56f5baSLois Curfman McInnes 324ea56f5baSLois Curfman McInnes Level: advanced 3258cbba510SBarry Smith 32629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 32729e0b56fSBarry Smith 32829e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 32929e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 33029e0b56fSBarry Smith @*/ 331af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 332329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 33329e0b56fSBarry Smith { 33429e0b56fSBarry Smith int ierr; 33529e0b56fSBarry Smith Scalar mone = -1.0; 3366831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3378f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 33829e0b56fSBarry Smith 3393a40ed3dSBarry Smith PetscFunctionBegin; 3408cbba510SBarry Smith *flag = 0; 34129e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 34229e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3438f99978dSLois Curfman McInnes if (neP->CheckStep) { 3448f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3458f99978dSLois Curfman McInnes } 346329f5518SBarry Smith 347329f5518SBarry Smith /* don't evaluate function the last time through */ 348329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 34929e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 350329f5518SBarry Smith } 35129e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3523a40ed3dSBarry Smith PetscFunctionReturn(0); 35329e0b56fSBarry Smith } 35404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 35529e0b56fSBarry Smith #undef __FUNC__ 356b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 3574b828684SBarry Smith /*@C 358f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3595e76c431SBarry Smith 360c7afd0dbSLois Curfman McInnes Collective on SNES 361c7afd0dbSLois Curfman McInnes 3625e76c431SBarry Smith Input Parameters: 363c7afd0dbSLois Curfman McInnes + snes - nonlinear context 364af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3655e76c431SBarry Smith . x - current iterate 3665e76c431SBarry Smith . f - residual evaluated at x 3675e76c431SBarry Smith . y - search direction (contains new iterate on output) 3685e76c431SBarry Smith . w - work vector 369c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3705e76c431SBarry Smith 371393d2d9aSLois Curfman McInnes Output Parameters: 372c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3735e76c431SBarry Smith . y - new iterate (contains search direction on input) 3745e76c431SBarry Smith . gnorm - 2-norm of g 3755e76c431SBarry Smith . ynorm - 2-norm of search length 376c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 377fee21e36SBarry Smith 378c4a48953SLois Curfman McInnes Options Database Key: 379c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 380c4a48953SLois Curfman McInnes 3815e76c431SBarry Smith Notes: 3825e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 3835e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 3845e76c431SBarry Smith 38536851e7fSLois Curfman McInnes Level: advanced 38636851e7fSLois Curfman McInnes 38728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 38828ae5a21SLois Curfman McInnes 389af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3905e76c431SBarry Smith @*/ 391af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 392329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3935e76c431SBarry Smith { 394406556e6SLois Curfman McInnes /* 395406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 396406556e6SLois Curfman McInnes minimization problem: 397406556e6SLois Curfman McInnes min z(x): R^n -> R, 398406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 399406556e6SLois Curfman McInnes */ 400406556e6SLois Curfman McInnes 401329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 402329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 403aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4045e42470aSBarry Smith Scalar cinitslope,clambda; 4056b5873e3SBarry Smith #endif 4065e42470aSBarry Smith int ierr,count; 4076831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4085334005bSBarry Smith Scalar mone = -1.0,scale; 4098f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4105e76c431SBarry Smith 4113a40ed3dSBarry Smith PetscFunctionBegin; 4127857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 413761aaf1bSLois Curfman McInnes *flag = 0; 4145e76c431SBarry Smith alpha = neP->alpha; 4155e76c431SBarry Smith maxstep = neP->maxstep; 4165e76c431SBarry Smith steptol = neP->steptol; 4175e76c431SBarry Smith 418cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 419a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 420a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 421a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 422a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 423a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 424ad922ac9SBarry Smith goto theend1; 42594a9d846SBarry Smith } 4265e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4275e42470aSBarry Smith scale = maxstep/(*ynorm); 428aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 429*90cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 4306b5873e3SBarry Smith #else 431*90cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 4326b5873e3SBarry Smith #endif 433393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4345e76c431SBarry Smith *ynorm = maxstep; 4355e76c431SBarry Smith } 4365e76c431SBarry Smith minlambda = steptol/(*ynorm); 437a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 438aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 439a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 440329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 4415e42470aSBarry Smith #else 442a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4435e42470aSBarry Smith #endif 4445e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4455e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4465e76c431SBarry Smith 447393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4485334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 44978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 450313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 451406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 452393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 45394a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 454ad922ac9SBarry Smith goto theend1; 4555e76c431SBarry Smith } 4565e76c431SBarry Smith 4575e76c431SBarry Smith /* Fit points with quadratic */ 458313b5538SBarry Smith lambda = 1.0; 459a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4605e76c431SBarry Smith lambdaprev = lambda; 4615e76c431SBarry Smith gnormprev = *gnorm; 462329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 463ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 464ddd12767SBarry Smith else lambda = lambdatemp; 465393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 466ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 467aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 468ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 4695e42470aSBarry Smith #else 470ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 4715e42470aSBarry Smith #endif 47278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 473cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 474406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 475393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 476ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 477ad922ac9SBarry Smith goto theend1; 4785e76c431SBarry Smith } 4795e76c431SBarry Smith 4805e76c431SBarry Smith /* Fit points with cubic */ 4815e76c431SBarry Smith count = 1; 4825e76c431SBarry Smith while (1) { 4835e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 4842b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 485*90cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 486393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 487761aaf1bSLois Curfman McInnes *flag = -1; break; 4885e76c431SBarry Smith } 489406556e6SLois Curfman McInnes t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 490406556e6SLois Curfman McInnes t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 491ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4922b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 4935e76c431SBarry Smith d = b*b - 3*a*initslope; 4945e76c431SBarry Smith if (d < 0.0) d = 0.0; 4955e76c431SBarry Smith if (a == 0.0) { 4965e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 4975e76c431SBarry Smith } else { 4985e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 4995e76c431SBarry Smith } 5005e76c431SBarry Smith lambdaprev = lambda; 5015e76c431SBarry Smith gnormprev = *gnorm; 502329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 503329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5045e76c431SBarry Smith else lambda = lambdatemp; 505393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 506ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 507aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 508ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 509393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5105e42470aSBarry Smith #else 511ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5125e42470aSBarry Smith #endif 51378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 514cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 515406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 516393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 517ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5182715565aSLois Curfman McInnes goto theend1; 5192b022350SLois Curfman McInnes } else { 5202b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5215e76c431SBarry Smith } 5225e76c431SBarry Smith count++; 5235e76c431SBarry Smith } 524ad922ac9SBarry Smith theend1: 5258f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5268f99978dSLois Curfman McInnes if (neP->CheckStep) { 5278f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5288f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5298f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5308f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5318f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5328f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5338f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5348f99978dSLois Curfman McInnes } 5358f99978dSLois Curfman McInnes } 5367857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5373a40ed3dSBarry Smith PetscFunctionReturn(0); 5385e76c431SBarry Smith } 53904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5405615d1e5SSatish Balay #undef __FUNC__ 541b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 5424b828684SBarry Smith /*@C 543f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5445e76c431SBarry Smith 545c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 546c7afd0dbSLois Curfman McInnes 5475e42470aSBarry Smith Input Parameters: 548c7afd0dbSLois Curfman McInnes + snes - the SNES context 549af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5505e42470aSBarry Smith . x - current iterate 5515e42470aSBarry Smith . f - residual evaluated at x 5525e42470aSBarry Smith . y - search direction (contains new iterate on output) 5535e42470aSBarry Smith . w - work vector 554c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5555e42470aSBarry Smith 556c4a48953SLois Curfman McInnes Output Parameters: 557c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5585e42470aSBarry Smith . y - new iterate (contains search direction on input) 5595e42470aSBarry Smith . gnorm - 2-norm of g 5605e42470aSBarry Smith . ynorm - 2-norm of search length 561c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 562fee21e36SBarry Smith 563c4a48953SLois Curfman McInnes Options Database Key: 564c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 565c4a48953SLois Curfman McInnes 5665e42470aSBarry Smith Notes: 5676831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 5685e42470aSBarry Smith 56936851e7fSLois Curfman McInnes Level: advanced 57036851e7fSLois Curfman McInnes 571f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 572f59ffdedSLois Curfman McInnes 573af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 5745e42470aSBarry Smith @*/ 575af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 576329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 5775e76c431SBarry Smith { 578406556e6SLois Curfman McInnes /* 579406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 580406556e6SLois Curfman McInnes minimization problem: 581406556e6SLois Curfman McInnes min z(x): R^n -> R, 582406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 583406556e6SLois Curfman McInnes */ 584329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 585aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5865e42470aSBarry Smith Scalar cinitslope,clambda; 5876b5873e3SBarry Smith #endif 5885e42470aSBarry Smith int ierr,count; 5896831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 5905334005bSBarry Smith Scalar mone = -1.0,scale; 5918f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 5925e76c431SBarry Smith 5933a40ed3dSBarry Smith PetscFunctionBegin; 5947857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 595761aaf1bSLois Curfman McInnes *flag = 0; 5965e42470aSBarry Smith alpha = neP->alpha; 5975e42470aSBarry Smith maxstep = neP->maxstep; 5985e42470aSBarry Smith steptol = neP->steptol; 5995e76c431SBarry Smith 6003505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 601b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 60294a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 603b37302c6SLois Curfman McInnes *gnorm = fnorm; 604b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 605b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 606ad922ac9SBarry Smith goto theend2; 60794a9d846SBarry Smith } 6085e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6095e42470aSBarry Smith scale = maxstep/(*ynorm); 610393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6115e42470aSBarry Smith *ynorm = maxstep; 6125e76c431SBarry Smith } 6135e42470aSBarry Smith minlambda = steptol/(*ynorm); 614a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 615aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 616a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 617329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6185e42470aSBarry Smith #else 619a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6205e42470aSBarry Smith #endif 6215e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6225e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6235e42470aSBarry Smith 624393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6255334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 62678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 627cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 628406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 629393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 63094a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 631ad922ac9SBarry Smith goto theend2; 6325e42470aSBarry Smith } 6335e42470aSBarry Smith 6345e42470aSBarry Smith /* Fit points with quadratic */ 635313b5538SBarry Smith lambda = 1.0; 6365e42470aSBarry Smith count = 1; 6375e42470aSBarry Smith while (1) { 6385e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 639981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 640981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 641ea4d3ed3SLois Curfman McInnes fnorm,*gnorm,*ynorm,lambda,initslope); 642393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 643761aaf1bSLois Curfman McInnes *flag = -1; break; 6445e42470aSBarry Smith } 645a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 646329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 647329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 648329f5518SBarry Smith else lambda = lambdatemp; 649393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6503505fcc1SBarry Smith lambdaneg = -lambda; 651aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6523505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6535e42470aSBarry Smith #else 6543505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6555e42470aSBarry Smith #endif 65678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 657cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 658406556e6SLois Curfman McInnes if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 659393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 660981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 66106259719SBarry Smith break; 6625e42470aSBarry Smith } 6635e42470aSBarry Smith count++; 6645e42470aSBarry Smith } 665ad922ac9SBarry Smith theend2: 6668f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6678f99978dSLois Curfman McInnes if (neP->CheckStep) { 6688f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6698f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6708f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6718f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6728f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6738f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6748f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6758f99978dSLois Curfman McInnes } 6768f99978dSLois Curfman McInnes } 6777857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6783a40ed3dSBarry Smith PetscFunctionReturn(0); 6795e76c431SBarry Smith } 68004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6815615d1e5SSatish Balay #undef __FUNC__ 682b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 683c9e6a524SLois Curfman McInnes /*@C 68477c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 6856831982aSBarry Smith by the method SNESEQLS. 6865e76c431SBarry Smith 6875e76c431SBarry Smith Input Parameters: 688c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 689af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 690c7afd0dbSLois Curfman McInnes - func - pointer to int function 6915e76c431SBarry Smith 692fee21e36SBarry Smith Collective on SNES 693fee21e36SBarry Smith 694c4a48953SLois Curfman McInnes Available Routines: 695c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 696f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 697f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 698af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 6995e76c431SBarry Smith 700c4a48953SLois Curfman McInnes Options Database Keys: 701af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 702c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 703c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 704c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 705c4a48953SLois Curfman McInnes 7065e76c431SBarry Smith Calling sequence of func: 707bd208895SLois Curfman McInnes .vb 708af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 709329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 710bd208895SLois Curfman McInnes .ve 7115e76c431SBarry Smith 7125e76c431SBarry Smith Input parameters for func: 713c7afd0dbSLois Curfman McInnes + snes - nonlinear context 714af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7155e76c431SBarry Smith . x - current iterate 7165e76c431SBarry Smith . f - residual evaluated at x 7175e76c431SBarry Smith . y - search direction (contains new iterate on output) 7185e76c431SBarry Smith . w - work vector 719c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7205e76c431SBarry Smith 7215e76c431SBarry Smith Output parameters for func: 722c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7235e76c431SBarry Smith . y - new iterate (contains search direction on input) 7245e76c431SBarry Smith . gnorm - 2-norm of g 7255e76c431SBarry Smith . ynorm - 2-norm of search length 726c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 727761aaf1bSLois Curfman McInnes on failure. 728f59ffdedSLois Curfman McInnes 72936851e7fSLois Curfman McInnes Level: advanced 73036851e7fSLois Curfman McInnes 731f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 732f59ffdedSLois Curfman McInnes 7334ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7345e76c431SBarry Smith @*/ 735329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7365e76c431SBarry Smith { 737329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 73882bf6240SBarry Smith 7393a40ed3dSBarry Smith PetscFunctionBegin; 74094d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 74182bf6240SBarry Smith if (f) { 742af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 74382bf6240SBarry Smith } 7443a40ed3dSBarry Smith PetscFunctionReturn(0); 7455e76c431SBarry Smith } 74604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 747fb2e594dSBarry Smith EXTERN_C_BEGIN 74882bf6240SBarry Smith #undef __FUNC__ 749b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 750af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 751af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 75282bf6240SBarry Smith { 75382bf6240SBarry Smith PetscFunctionBegin; 7546831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7556831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 75682bf6240SBarry Smith PetscFunctionReturn(0); 75782bf6240SBarry Smith } 758fb2e594dSBarry Smith EXTERN_C_END 75904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 760c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 761b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 762c8dd0c0dSLois Curfman McInnes /*@C 763530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7646831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 765c8dd0c0dSLois Curfman McInnes 766c8dd0c0dSLois Curfman McInnes Input Parameters: 767c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 768c8dd0c0dSLois Curfman McInnes . func - pointer to int function 769c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 770c8dd0c0dSLois Curfman McInnes 771c8dd0c0dSLois Curfman McInnes Collective on SNES 772c8dd0c0dSLois Curfman McInnes 773c8dd0c0dSLois Curfman McInnes Calling sequence of func: 774c8dd0c0dSLois Curfman McInnes .vb 775b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 776c8dd0c0dSLois Curfman McInnes .ve 777b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 778b1ae27eaSLois Curfman McInnes on failure. 779c8dd0c0dSLois Curfman McInnes 780c8dd0c0dSLois Curfman McInnes Input parameters for func: 781c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 782c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 7838f99978dSLois Curfman McInnes - x - current candidate iterate 784c8dd0c0dSLois Curfman McInnes 785c8dd0c0dSLois Curfman McInnes Output parameters for func: 786c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 787c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 788c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 789c8dd0c0dSLois Curfman McInnes 790c8dd0c0dSLois Curfman McInnes Level: advanced 791c8dd0c0dSLois Curfman McInnes 7928f99978dSLois Curfman McInnes Notes: 793b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 794b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 795b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 796b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 797ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 7988f99978dSLois Curfman McInnes 799b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 800b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 801b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 802ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 803ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 804ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 805ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 806b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8078f99978dSLois Curfman McInnes 808c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 809c8dd0c0dSLois Curfman McInnes 810c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 811c8dd0c0dSLois Curfman McInnes @*/ 8128f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 813c8dd0c0dSLois Curfman McInnes { 8148f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 815c8dd0c0dSLois Curfman McInnes 816c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 817c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 818c8dd0c0dSLois Curfman McInnes if (f) { 819c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 820c8dd0c0dSLois Curfman McInnes } 821c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 822c8dd0c0dSLois Curfman McInnes } 823c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 824c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 825c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 826b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8278f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 828c8dd0c0dSLois Curfman McInnes { 829c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8306831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8316831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 832c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 833c8dd0c0dSLois Curfman McInnes } 834c8dd0c0dSLois Curfman McInnes EXTERN_C_END 835c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 83604d965bbSLois Curfman McInnes /* 8376831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 83882bf6240SBarry Smith 83904d965bbSLois Curfman McInnes Input Parameter: 84004d965bbSLois Curfman McInnes . snes - the SNES context 84104d965bbSLois Curfman McInnes 84204d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 84304d965bbSLois Curfman McInnes */ 8445615d1e5SSatish Balay #undef __FUNC__ 845b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 846f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8475e42470aSBarry Smith { 8486831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 849d132466eSBarry Smith int ierr; 8506daaf66cSBarry Smith 8513a40ed3dSBarry Smith PetscFunctionBegin; 8526831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 853d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 854d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 855d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 856d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8573a40ed3dSBarry Smith PetscFunctionReturn(0); 8585e42470aSBarry Smith } 85904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 86004d965bbSLois Curfman McInnes /* 8616831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 86204d965bbSLois Curfman McInnes 86304d965bbSLois Curfman McInnes Input Parameters: 86404d965bbSLois Curfman McInnes . SNES - the SNES context 86504d965bbSLois Curfman McInnes . viewer - visualization context 86604d965bbSLois Curfman McInnes 86704d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 86804d965bbSLois Curfman McInnes */ 8695615d1e5SSatish Balay #undef __FUNC__ 870b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 871e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 872a935fc98SLois Curfman McInnes { 8736831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 87419bcc07fSBarry Smith char *cstr; 87551695f54SLois Curfman McInnes int ierr; 8766831982aSBarry Smith PetscTruth isascii; 877a935fc98SLois Curfman McInnes 8783a40ed3dSBarry Smith PetscFunctionBegin; 8796831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8800f5bd95cSBarry Smith if (isascii) { 88119bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 88219bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 88319bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 88419bcc07fSBarry Smith else cstr = "unknown"; 885d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 886d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 8875cd90555SBarry Smith } else { 8880f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 88919bcc07fSBarry Smith } 8903a40ed3dSBarry Smith PetscFunctionReturn(0); 891a935fc98SLois Curfman McInnes } 89204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 89304d965bbSLois Curfman McInnes /* 8946831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 89504d965bbSLois Curfman McInnes 89604d965bbSLois Curfman McInnes Input Parameter: 89704d965bbSLois Curfman McInnes . snes - the SNES context 89804d965bbSLois Curfman McInnes 89904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 90004d965bbSLois Curfman McInnes */ 9015615d1e5SSatish Balay #undef __FUNC__ 902b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 903f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 9045e42470aSBarry Smith { 9056831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 9065e42470aSBarry Smith char ver[16]; 9075e42470aSBarry Smith double tmp; 908f1af5d2fSBarry Smith int ierr; 909f1af5d2fSBarry Smith PetscTruth flg; 9105e42470aSBarry Smith 9113a40ed3dSBarry Smith PetscFunctionBegin; 91209d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 91325cdf11fSBarry Smith if (flg) { 9145e42470aSBarry Smith ls->alpha = tmp; 9155e42470aSBarry Smith } 91609d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 91725cdf11fSBarry Smith if (flg) { 9185e42470aSBarry Smith ls->maxstep = tmp; 9195e42470aSBarry Smith } 92009d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 92125cdf11fSBarry Smith if (flg) { 9225e42470aSBarry Smith ls->steptol = tmp; 9235e42470aSBarry Smith } 92409d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 92525cdf11fSBarry Smith if (flg) { 926c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9270f5bd95cSBarry Smith 928c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 929c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 930c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 931c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9320f5bd95cSBarry Smith 9330f5bd95cSBarry Smith if (isbasic) { 934af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9350f5bd95cSBarry Smith } else if (isnonorms) { 936af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9370f5bd95cSBarry Smith } else if (isquad) { 938af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9390f5bd95cSBarry Smith } else if (iscubic) { 940af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9415e42470aSBarry Smith } 942a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9435e42470aSBarry Smith } 9443a40ed3dSBarry Smith PetscFunctionReturn(0); 9455e42470aSBarry Smith } 94604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 94704d965bbSLois Curfman McInnes /* 9486831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9496831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 95004d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 95104d965bbSLois Curfman McInnes 95204d965bbSLois Curfman McInnes 95304d965bbSLois Curfman McInnes Input Parameter: 95404d965bbSLois Curfman McInnes . snes - the SNES context 95504d965bbSLois Curfman McInnes 95604d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 95704d965bbSLois Curfman McInnes */ 958fb2e594dSBarry Smith EXTERN_C_BEGIN 9595615d1e5SSatish Balay #undef __FUNC__ 960b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 961f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9625e42470aSBarry Smith { 96382bf6240SBarry Smith int ierr; 9646831982aSBarry Smith SNES_EQ_LS *neP; 9655e42470aSBarry Smith 9663a40ed3dSBarry Smith PetscFunctionBegin; 967a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 968a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 969a8c6a408SBarry Smith } 97082bf6240SBarry Smith 971f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 972f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 973f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 974f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 975f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 976f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 977f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 9785baf8537SBarry Smith snes->nwork = 0; 9795e42470aSBarry Smith 9806831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 9816831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 9825e42470aSBarry Smith snes->data = (void*)neP; 9835e42470aSBarry Smith neP->alpha = 1.e-4; 9845e42470aSBarry Smith neP->maxstep = 1.e8; 9855e42470aSBarry Smith neP->steptol = 1.e-12; 9865e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 987c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 988c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 989c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 99082bf6240SBarry Smith 991f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 9920c97f09dSSatish Balay SNESSetLineSearch_LS);CHKERRQ(ierr); 993f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 9940c97f09dSSatish Balay SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 99582bf6240SBarry Smith 9963a40ed3dSBarry Smith PetscFunctionReturn(0); 9975e42470aSBarry Smith } 998fb2e594dSBarry Smith EXTERN_C_END 99982bf6240SBarry Smith 100082bf6240SBarry Smith 100182bf6240SBarry Smith 1002