1*85471664SBarry Smith /*$Id: ls.c,v 1.160 2000/08/01 20:57:22 bsmith Exp bsmith $*/ 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 58a5d9ee4SBarry Smith /* 68a5d9ee4SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the function, 78a5d9ee4SBarry Smith but not a zero. In the case when one cannot compute J^T F we use the fact that 836669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 936669109SBarry Smith for this trick. 108a5d9ee4SBarry Smith */ 118a5d9ee4SBarry Smith #undef __FUNC__ 128a5d9ee4SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"SNESLSCheckLocalMin_Private" 138a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 168a5d9ee4SBarry Smith int ierr; 1736669109SBarry Smith PetscTruth hastranspose; 188a5d9ee4SBarry Smith 198a5d9ee4SBarry Smith PetscFunctionBegin; 208a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2136669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2236669109SBarry Smith if (hastranspose) { 238a5d9ee4SBarry Smith /* Compute || J^T F|| */ 248a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 258a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 265d1a10b1SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 3036669109SBarry Smith Scalar result; 3136669109SBarry Smith PetscReal wnorm; 3236669109SBarry Smith 3336669109SBarry Smith ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDestroy(work);CHKERRQ(ierr); 3936669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 405d1a10b1SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 4774637425SBarry Smith Checks if J^T(F - AX) = 0 4874637425SBarry Smith */ 4974637425SBarry Smith #undef __FUNC__ 5074637425SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckResidual_Private"></a>*/"SNESLSCheckResidual_Private" 5174637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5274637425SBarry Smith { 5374637425SBarry Smith PetscReal a1,a2; 5474637425SBarry Smith int ierr; 5574637425SBarry Smith PetscTruth hastranspose; 5674637425SBarry Smith Scalar mone = -1.0; 5774637425SBarry Smith 5874637425SBarry Smith PetscFunctionBegin; 5974637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 6074637425SBarry Smith if (hastranspose) { 6174637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6274637425SBarry Smith ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr); 6374637425SBarry Smith 6474637425SBarry Smith /* Compute || J^T W|| */ 6574637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6774637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 68*85471664SBarry Smith if (PetscAbsScalar(a1) != 0) { 69*85471664SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1); 70*85471664SBarry Smith } 7174637425SBarry Smith } 7274637425SBarry Smith PetscFunctionReturn(0); 7374637425SBarry Smith } 7474637425SBarry Smith 7504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7604d965bbSLois Curfman McInnes 7704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7804d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 7904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8004d965bbSLois Curfman McInnes respectively. 8104d965bbSLois Curfman McInnes 82fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 85fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 8804d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 8904d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 9004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9304d965bbSLois Curfman McInnes for all nonlinear solvers. 9404d965bbSLois Curfman McInnes 9504d965bbSLois Curfman McInnes Another key routine is: 9604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9904d965bbSLois Curfman McInnes SNESSolve() if necessary. 10004d965bbSLois Curfman McInnes 10104d965bbSLois Curfman McInnes Additional basic routines are: 10204d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 10304d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10404d965bbSLois Curfman McInnes have actually been used. 10504d965bbSLois Curfman McInnes These are called by application codes via the interface routines 10604d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 10704d965bbSLois Curfman McInnes 10804d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10904d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 11004d965bbSLois Curfman McInnes above description applies to these categories also. 11104d965bbSLois Curfman McInnes 11204d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1135e42470aSBarry Smith /* 11404d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 11504d965bbSLois Curfman McInnes method with a line search. 1165e76c431SBarry Smith 11704d965bbSLois Curfman McInnes Input Parameters: 11804d965bbSLois Curfman McInnes . snes - the SNES context 1195e76c431SBarry Smith 12004d965bbSLois Curfman McInnes Output Parameter: 121c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12204d965bbSLois Curfman McInnes 12304d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1245e76c431SBarry Smith 1255e76c431SBarry Smith Notes: 1265e76c431SBarry Smith This implements essentially a truncated Newton method with a 1275e76c431SBarry Smith line search. By default a cubic backtracking line search 1285e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1295e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 130393d2d9aSLois Curfman McInnes and Schnabel. 1315e76c431SBarry Smith */ 1325615d1e5SSatish Balay #undef __FUNC__ 133b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 134f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 1355e76c431SBarry Smith { 1366831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 13784c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 138112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1405e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 1415e76c431SBarry Smith 1423a40ed3dSBarry Smith PetscFunctionBegin; 143184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 144184914b5SBarry Smith 1455e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1465e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14739e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1485e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1495e42470aSBarry Smith G = snes->work[1]; 1505e42470aSBarry Smith W = snes->work[2]; 1515e76c431SBarry Smith 1524c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1534c49b128SBarry Smith snes->iter = 0; 1544c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1555334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 156cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1570f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1585e42470aSBarry Smith snes->norm = fnorm; 1590f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16084c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16194a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1625e76c431SBarry Smith 163184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16494a9d846SBarry Smith 165d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 166d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 167d034289bSBarry Smith 1685e76c431SBarry Smith for (i=0; i<maxits; i++) { 1695e76c431SBarry Smith 170ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1715334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1725334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 17378b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 17474637425SBarry Smith 175*85471664SBarry Smith if (PLogPrintInfo){ 17674637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 177*85471664SBarry Smith } 17874637425SBarry Smith 17990cbd584SBarry Smith /* should check what happened to the linear solve? */ 1803505fcc1SBarry Smith snes->linear_its += lits; 181af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 182ea4d3ed3SLois Curfman McInnes 183ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 184ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 185ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 186ea4d3ed3SLois Curfman McInnes */ 18781b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 188af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 189af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 1903505fcc1SBarry Smith if (lsfail) { 1918a5d9ee4SBarry Smith PetscTruth ismin; 1923505fcc1SBarry Smith snes->nfailures++; 1933505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 1948a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 1958a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1963505fcc1SBarry Smith break; 1973505fcc1SBarry Smith } 1985e76c431SBarry Smith 19939e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 20039e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 2015e76c431SBarry Smith fnorm = gnorm; 2025e76c431SBarry Smith 2030f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 20484c569c9SBarry Smith snes->iter = i+1; 2055e42470aSBarry Smith snes->norm = fnorm; 2060f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 20784c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 20894a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 2095e76c431SBarry Smith 2105e76c431SBarry Smith /* Test for convergence */ 21129e0b56fSBarry Smith if (snes->converged) { 21229e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2133505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2143505fcc1SBarry Smith if (snes->reason) { 21516c95cacSBarry Smith break; 21616c95cacSBarry Smith } 21716c95cacSBarry Smith } 21829e0b56fSBarry Smith } 21939e2f89bSBarry Smith if (X != snes->vec_sol) { 220393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 221e15f7bb5SBarry Smith } 222e15f7bb5SBarry Smith if (F != snes->vec_func) { 223e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 224e15f7bb5SBarry Smith } 22539e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 22639e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 22752392280SLois Curfman McInnes if (i == maxits) { 228981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 22952392280SLois Curfman McInnes i--; 2303505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 23152392280SLois Curfman McInnes } 2320f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2330f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2345e42470aSBarry Smith *outits = i+1; 2353a40ed3dSBarry Smith PetscFunctionReturn(0); 2365e76c431SBarry Smith } 23704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 23804d965bbSLois Curfman McInnes /* 23904d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 2406831982aSBarry Smith of the SNESEQLS nonlinear solver. 24104d965bbSLois Curfman McInnes 24204d965bbSLois Curfman McInnes Input Parameter: 24304d965bbSLois Curfman McInnes . snes - the SNES context 24404d965bbSLois Curfman McInnes . x - the solution vector 24504d965bbSLois Curfman McInnes 24604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 24704d965bbSLois Curfman McInnes 24804d965bbSLois Curfman McInnes Notes: 24904d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 25004d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 25104d965bbSLois Curfman McInnes the call to SNESSolve(). 25204d965bbSLois Curfman McInnes */ 2535615d1e5SSatish Balay #undef __FUNC__ 254b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 255f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 2565e76c431SBarry Smith { 2575e42470aSBarry Smith int ierr; 2583a40ed3dSBarry Smith 2593a40ed3dSBarry Smith PetscFunctionBegin; 26081b6cf68SLois Curfman McInnes snes->nwork = 4; 261d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 2625e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 26381b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2643a40ed3dSBarry Smith PetscFunctionReturn(0); 2655e76c431SBarry Smith } 26604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 26704d965bbSLois Curfman McInnes /* 2686831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 26904d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 27004d965bbSLois Curfman McInnes 27104d965bbSLois Curfman McInnes Input Parameter: 27204d965bbSLois Curfman McInnes . snes - the SNES context 27304d965bbSLois Curfman McInnes 274de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 27504d965bbSLois Curfman McInnes */ 2765615d1e5SSatish Balay #undef __FUNC__ 277b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 278e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2795e76c431SBarry Smith { 280393d2d9aSLois Curfman McInnes int ierr; 2813a40ed3dSBarry Smith 2823a40ed3dSBarry Smith PetscFunctionBegin; 2835baf8537SBarry Smith if (snes->nwork) { 2844b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 28521c89e3eSBarry Smith } 2865c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2873a40ed3dSBarry Smith PetscFunctionReturn(0); 2885e76c431SBarry Smith } 28904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2905615d1e5SSatish Balay #undef __FUNC__ 291b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 29282bf6240SBarry Smith 2934b828684SBarry Smith /*@C 2945e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2955e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2965e76c431SBarry Smith to serve as a template and is not recommended for general use. 2975e76c431SBarry Smith 298c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 299c7afd0dbSLois Curfman McInnes 3005e76c431SBarry Smith Input Parameters: 301c7afd0dbSLois Curfman McInnes + snes - nonlinear context 302af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3035e76c431SBarry Smith . x - current iterate 3045e76c431SBarry Smith . f - residual evaluated at x 3055e76c431SBarry Smith . y - search direction (contains new iterate on output) 3065e76c431SBarry Smith . w - work vector 307c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3085e76c431SBarry Smith 309c4a48953SLois Curfman McInnes Output Parameters: 310c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3115e76c431SBarry Smith . y - new iterate (contains search direction on input) 3125e76c431SBarry Smith . gnorm - 2-norm of g 3135e76c431SBarry Smith . ynorm - 2-norm of search length 314c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 315fee21e36SBarry Smith 316c4a48953SLois Curfman McInnes Options Database Key: 317c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 318c4a48953SLois Curfman McInnes 31936851e7fSLois Curfman McInnes Level: advanced 32036851e7fSLois Curfman McInnes 32128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 32228ae5a21SLois Curfman McInnes 323f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 324af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3255e76c431SBarry Smith @*/ 3265d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3275e76c431SBarry Smith { 3285e42470aSBarry Smith int ierr; 3295334005bSBarry Smith Scalar mone = -1.0; 3306831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3318f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3325334005bSBarry Smith 3333a40ed3dSBarry Smith PetscFunctionBegin; 334761aaf1bSLois Curfman McInnes *flag = 0; 3357857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 336a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 337ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3388f99978dSLois Curfman McInnes if (neP->CheckStep) { 3398f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3408f99978dSLois Curfman McInnes } 341ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 342a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 3437857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3443a40ed3dSBarry Smith PetscFunctionReturn(0); 3455e76c431SBarry Smith } 34604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3475615d1e5SSatish Balay #undef __FUNC__ 348b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 34982bf6240SBarry Smith 35029e0b56fSBarry Smith /*@C 35129e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 35229e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 35329e0b56fSBarry Smith even compute the norm of the function or search direction; this 35429e0b56fSBarry Smith is intended only when you know the full step is fine and are 35529e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 356c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 357c7afd0dbSLois Curfman McInnes 358c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 35929e0b56fSBarry Smith 36029e0b56fSBarry Smith Input Parameters: 361c7afd0dbSLois Curfman McInnes + snes - nonlinear context 362af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 36329e0b56fSBarry Smith . x - current iterate 36429e0b56fSBarry Smith . f - residual evaluated at x 36529e0b56fSBarry Smith . y - search direction (contains new iterate on output) 36629e0b56fSBarry Smith . w - work vector 367c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 36829e0b56fSBarry Smith 36929e0b56fSBarry Smith Output Parameters: 370c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 37129e0b56fSBarry Smith . gnorm - not changed 37229e0b56fSBarry Smith . ynorm - not changed 373c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 374fee21e36SBarry Smith 37529e0b56fSBarry Smith Options Database Key: 376c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 37729e0b56fSBarry Smith 3788cbba510SBarry Smith Notes: 379ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 380ea56f5baSLois Curfman McInnes either the options 381ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 382ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 383329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 384329f5518SBarry Smith otherwise, the SNES solver will generate an error. 385329f5518SBarry Smith 386329f5518SBarry Smith During the final iteration this will not evaluate the function at 387329f5518SBarry Smith the solution point. This is to save a function evaluation while 388329f5518SBarry Smith using pseudo-timestepping. 3898cbba510SBarry Smith 390ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 391ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 392ea56f5baSLois Curfman McInnes correct, since they are not computed. 393ea56f5baSLois Curfman McInnes 394ea56f5baSLois Curfman McInnes Level: advanced 3958cbba510SBarry Smith 39629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 39729e0b56fSBarry Smith 39829e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 39929e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 40029e0b56fSBarry Smith @*/ 4015d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 40229e0b56fSBarry Smith { 40329e0b56fSBarry Smith int ierr; 40429e0b56fSBarry Smith Scalar mone = -1.0; 4056831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4068f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 40729e0b56fSBarry Smith 4083a40ed3dSBarry Smith PetscFunctionBegin; 4098cbba510SBarry Smith *flag = 0; 41029e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 41129e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4128f99978dSLois Curfman McInnes if (neP->CheckStep) { 4138f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4148f99978dSLois Curfman McInnes } 415329f5518SBarry Smith 416329f5518SBarry Smith /* don't evaluate function the last time through */ 417329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 41829e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 419329f5518SBarry Smith } 42029e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4213a40ed3dSBarry Smith PetscFunctionReturn(0); 42229e0b56fSBarry Smith } 42304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 42429e0b56fSBarry Smith #undef __FUNC__ 425b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 4264b828684SBarry Smith /*@C 427f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4285e76c431SBarry Smith 429c7afd0dbSLois Curfman McInnes Collective on SNES 430c7afd0dbSLois Curfman McInnes 4315e76c431SBarry Smith Input Parameters: 432c7afd0dbSLois Curfman McInnes + snes - nonlinear context 433af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4345e76c431SBarry Smith . x - current iterate 4355e76c431SBarry Smith . f - residual evaluated at x 4365e76c431SBarry Smith . y - search direction (contains new iterate on output) 4375e76c431SBarry Smith . w - work vector 438c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4395e76c431SBarry Smith 440393d2d9aSLois Curfman McInnes Output Parameters: 441c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4425e76c431SBarry Smith . y - new iterate (contains search direction on input) 4435e76c431SBarry Smith . gnorm - 2-norm of g 4445e76c431SBarry Smith . ynorm - 2-norm of search length 445c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 446fee21e36SBarry Smith 447c4a48953SLois Curfman McInnes Options Database Key: 448c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 449c4a48953SLois Curfman McInnes 4505e76c431SBarry Smith Notes: 4515e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4525e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4535e76c431SBarry Smith 45436851e7fSLois Curfman McInnes Level: advanced 45536851e7fSLois Curfman McInnes 45628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 45728ae5a21SLois Curfman McInnes 458af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4595e76c431SBarry Smith @*/ 4605d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 4615e76c431SBarry Smith { 462406556e6SLois Curfman McInnes /* 463406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 464406556e6SLois Curfman McInnes minimization problem: 465406556e6SLois Curfman McInnes min z(x): R^n -> R, 466406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 467406556e6SLois Curfman McInnes */ 468406556e6SLois Curfman McInnes 469329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 470329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 471aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4725e42470aSBarry Smith Scalar cinitslope,clambda; 4736b5873e3SBarry Smith #endif 4745e42470aSBarry Smith int ierr,count; 4756831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4765334005bSBarry Smith Scalar mone = -1.0,scale; 4778f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4785e76c431SBarry Smith 4793a40ed3dSBarry Smith PetscFunctionBegin; 4807857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 481761aaf1bSLois Curfman McInnes *flag = 0; 4825e76c431SBarry Smith alpha = neP->alpha; 4835e76c431SBarry Smith maxstep = neP->maxstep; 4845e76c431SBarry Smith steptol = neP->steptol; 4855e76c431SBarry Smith 486cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 487a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 488a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 489a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 490a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 491a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 492ad922ac9SBarry Smith goto theend1; 49394a9d846SBarry Smith } 4945e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4955e42470aSBarry Smith scale = maxstep/(*ynorm); 496aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 49790cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 4986b5873e3SBarry Smith #else 49990cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 5006b5873e3SBarry Smith #endif 501393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5025e76c431SBarry Smith *ynorm = maxstep; 5035e76c431SBarry Smith } 5045e76c431SBarry Smith minlambda = steptol/(*ynorm); 505a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 506aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 507a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 508329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5095e42470aSBarry Smith #else 510a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5115e42470aSBarry Smith #endif 5125e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5135e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5145e76c431SBarry Smith 515393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5165334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 51778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 518313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5195d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 520393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 52194a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 522ad922ac9SBarry Smith goto theend1; 5235e76c431SBarry Smith } 5245e76c431SBarry Smith 5255e76c431SBarry Smith /* Fit points with quadratic */ 526313b5538SBarry Smith lambda = 1.0; 527a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5285e76c431SBarry Smith lambdaprev = lambda; 5295e76c431SBarry Smith gnormprev = *gnorm; 530329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 531ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 532ddd12767SBarry Smith else lambda = lambdatemp; 533393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 534ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 535aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 536ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5375e42470aSBarry Smith #else 538ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5395e42470aSBarry Smith #endif 54078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 541cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5425d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 543393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 544ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 545ad922ac9SBarry Smith goto theend1; 5465e76c431SBarry Smith } 5475e76c431SBarry Smith 5485e76c431SBarry Smith /* Fit points with cubic */ 5495e76c431SBarry Smith count = 1; 5505e76c431SBarry Smith while (1) { 5515e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 5522b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 55390cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 554393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 555761aaf1bSLois Curfman McInnes *flag = -1; break; 5565e76c431SBarry Smith } 5575d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 5585d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 559ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5602b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5615e76c431SBarry Smith d = b*b - 3*a*initslope; 5625e76c431SBarry Smith if (d < 0.0) d = 0.0; 5635e76c431SBarry Smith if (a == 0.0) { 5645e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5655e76c431SBarry Smith } else { 5665e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5675e76c431SBarry Smith } 5685e76c431SBarry Smith lambdaprev = lambda; 5695e76c431SBarry Smith gnormprev = *gnorm; 570329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 571329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5725e76c431SBarry Smith else lambda = lambdatemp; 573393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 574ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 575aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 576ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 577393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5785e42470aSBarry Smith #else 579ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5805e42470aSBarry Smith #endif 58178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 582cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5835d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 584393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 585ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5862715565aSLois Curfman McInnes goto theend1; 5872b022350SLois Curfman McInnes } else { 5882b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5895e76c431SBarry Smith } 5905e76c431SBarry Smith count++; 5915e76c431SBarry Smith } 592ad922ac9SBarry Smith theend1: 5938f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5948f99978dSLois Curfman McInnes if (neP->CheckStep) { 5958f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5968f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5978f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5988f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5998f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6008f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6018f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6028f99978dSLois Curfman McInnes } 6038f99978dSLois Curfman McInnes } 6047857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6053a40ed3dSBarry Smith PetscFunctionReturn(0); 6065e76c431SBarry Smith } 60704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6085615d1e5SSatish Balay #undef __FUNC__ 609b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 6104b828684SBarry Smith /*@C 611f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6125e76c431SBarry Smith 613c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 614c7afd0dbSLois Curfman McInnes 6155e42470aSBarry Smith Input Parameters: 616c7afd0dbSLois Curfman McInnes + snes - the SNES context 617af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6185e42470aSBarry Smith . x - current iterate 6195e42470aSBarry Smith . f - residual evaluated at x 6205e42470aSBarry Smith . y - search direction (contains new iterate on output) 6215e42470aSBarry Smith . w - work vector 622c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6235e42470aSBarry Smith 624c4a48953SLois Curfman McInnes Output Parameters: 625c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6265e42470aSBarry Smith . y - new iterate (contains search direction on input) 6275e42470aSBarry Smith . gnorm - 2-norm of g 6285e42470aSBarry Smith . ynorm - 2-norm of search length 629c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 630fee21e36SBarry Smith 631c4a48953SLois Curfman McInnes Options Database Key: 632c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 633c4a48953SLois Curfman McInnes 6345e42470aSBarry Smith Notes: 6356831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 6365e42470aSBarry Smith 63736851e7fSLois Curfman McInnes Level: advanced 63836851e7fSLois Curfman McInnes 639f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 640f59ffdedSLois Curfman McInnes 641af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6425e42470aSBarry Smith @*/ 6435d1a10b1SBarry Smith int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 6445e76c431SBarry Smith { 645406556e6SLois Curfman McInnes /* 646406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 647406556e6SLois Curfman McInnes minimization problem: 648406556e6SLois Curfman McInnes min z(x): R^n -> R, 649406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 650406556e6SLois Curfman McInnes */ 651329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 652aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6535e42470aSBarry Smith Scalar cinitslope,clambda; 6546b5873e3SBarry Smith #endif 6555e42470aSBarry Smith int ierr,count; 6566831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 6575334005bSBarry Smith Scalar mone = -1.0,scale; 6588f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6595e76c431SBarry Smith 6603a40ed3dSBarry Smith PetscFunctionBegin; 6617857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 662761aaf1bSLois Curfman McInnes *flag = 0; 6635e42470aSBarry Smith alpha = neP->alpha; 6645e42470aSBarry Smith maxstep = neP->maxstep; 6655e42470aSBarry Smith steptol = neP->steptol; 6665e76c431SBarry Smith 6673505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 668b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 66994a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 670b37302c6SLois Curfman McInnes *gnorm = fnorm; 671b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 672b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 673ad922ac9SBarry Smith goto theend2; 67494a9d846SBarry Smith } 6755e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6765e42470aSBarry Smith scale = maxstep/(*ynorm); 677393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6785e42470aSBarry Smith *ynorm = maxstep; 6795e76c431SBarry Smith } 6805e42470aSBarry Smith minlambda = steptol/(*ynorm); 681a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 682aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 683a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 684329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6855e42470aSBarry Smith #else 686a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6875e42470aSBarry Smith #endif 6885e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6895e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6905e42470aSBarry Smith 691393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6925334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 69378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 694cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 6955d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 696393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 69794a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 698ad922ac9SBarry Smith goto theend2; 6995e42470aSBarry Smith } 7005e42470aSBarry Smith 7015e42470aSBarry Smith /* Fit points with quadratic */ 702313b5538SBarry Smith lambda = 1.0; 7035e42470aSBarry Smith count = 1; 7045e42470aSBarry Smith while (1) { 7055e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 706981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 7075d1a10b1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 708393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 709761aaf1bSLois Curfman McInnes *flag = -1; break; 7105e42470aSBarry Smith } 711a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 712329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 713329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 714329f5518SBarry Smith else lambda = lambdatemp; 715393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7163505fcc1SBarry Smith lambdaneg = -lambda; 717aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7183505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7195e42470aSBarry Smith #else 7203505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7215e42470aSBarry Smith #endif 72278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 723cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7245d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 725393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 726981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 72706259719SBarry Smith break; 7285e42470aSBarry Smith } 7295e42470aSBarry Smith count++; 7305e42470aSBarry Smith } 731ad922ac9SBarry Smith theend2: 7328f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7338f99978dSLois Curfman McInnes if (neP->CheckStep) { 7348f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7358f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7368f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7378f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7388f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7398f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7408f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7418f99978dSLois Curfman McInnes } 7428f99978dSLois Curfman McInnes } 7437857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 7443a40ed3dSBarry Smith PetscFunctionReturn(0); 7455e76c431SBarry Smith } 74604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7475615d1e5SSatish Balay #undef __FUNC__ 748b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 749c9e6a524SLois Curfman McInnes /*@C 75077c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7516831982aSBarry Smith by the method SNESEQLS. 7525e76c431SBarry Smith 7535e76c431SBarry Smith Input Parameters: 754c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 755af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 756c7afd0dbSLois Curfman McInnes - func - pointer to int function 7575e76c431SBarry Smith 758fee21e36SBarry Smith Collective on SNES 759fee21e36SBarry Smith 760c4a48953SLois Curfman McInnes Available Routines: 761c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 762f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 763f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 764af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7655e76c431SBarry Smith 766c4a48953SLois Curfman McInnes Options Database Keys: 767af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 768c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 769c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 770c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 771c4a48953SLois Curfman McInnes 7725e76c431SBarry Smith Calling sequence of func: 773bd208895SLois Curfman McInnes .vb 774af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 775329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 776bd208895SLois Curfman McInnes .ve 7775e76c431SBarry Smith 7785e76c431SBarry Smith Input parameters for func: 779c7afd0dbSLois Curfman McInnes + snes - nonlinear context 780af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7815e76c431SBarry Smith . x - current iterate 7825e76c431SBarry Smith . f - residual evaluated at x 7835e76c431SBarry Smith . y - search direction (contains new iterate on output) 7845e76c431SBarry Smith . w - work vector 785c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7865e76c431SBarry Smith 7875e76c431SBarry Smith Output parameters for func: 788c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7895e76c431SBarry Smith . y - new iterate (contains search direction on input) 7905e76c431SBarry Smith . gnorm - 2-norm of g 7915e76c431SBarry Smith . ynorm - 2-norm of search length 792c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 793761aaf1bSLois Curfman McInnes on failure. 794f59ffdedSLois Curfman McInnes 79536851e7fSLois Curfman McInnes Level: advanced 79636851e7fSLois Curfman McInnes 797f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 798f59ffdedSLois Curfman McInnes 7995d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 8005d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 8015e76c431SBarry Smith @*/ 802329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 8035e76c431SBarry Smith { 804329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 80582bf6240SBarry Smith 8063a40ed3dSBarry Smith PetscFunctionBegin; 80794d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 80882bf6240SBarry Smith if (f) { 809af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 81082bf6240SBarry Smith } 8113a40ed3dSBarry Smith PetscFunctionReturn(0); 8125e76c431SBarry Smith } 81304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 814fb2e594dSBarry Smith EXTERN_C_BEGIN 81582bf6240SBarry Smith #undef __FUNC__ 816b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 817af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 818af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 81982bf6240SBarry Smith { 82082bf6240SBarry Smith PetscFunctionBegin; 8216831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 8226831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 82382bf6240SBarry Smith PetscFunctionReturn(0); 82482bf6240SBarry Smith } 825fb2e594dSBarry Smith EXTERN_C_END 82604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 827c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 828b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 829c8dd0c0dSLois Curfman McInnes /*@C 830530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 8316831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 832c8dd0c0dSLois Curfman McInnes 833c8dd0c0dSLois Curfman McInnes Input Parameters: 834c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 835c8dd0c0dSLois Curfman McInnes . func - pointer to int function 836c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 837c8dd0c0dSLois Curfman McInnes 838c8dd0c0dSLois Curfman McInnes Collective on SNES 839c8dd0c0dSLois Curfman McInnes 840c8dd0c0dSLois Curfman McInnes Calling sequence of func: 841c8dd0c0dSLois Curfman McInnes .vb 842b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 843c8dd0c0dSLois Curfman McInnes .ve 844b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 845b1ae27eaSLois Curfman McInnes on failure. 846c8dd0c0dSLois Curfman McInnes 847c8dd0c0dSLois Curfman McInnes Input parameters for func: 848c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 849c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8508f99978dSLois Curfman McInnes - x - current candidate iterate 851c8dd0c0dSLois Curfman McInnes 852c8dd0c0dSLois Curfman McInnes Output parameters for func: 853c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 854c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 855c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 856c8dd0c0dSLois Curfman McInnes 857c8dd0c0dSLois Curfman McInnes Level: advanced 858c8dd0c0dSLois Curfman McInnes 8598f99978dSLois Curfman McInnes Notes: 860b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 861b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 862b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 863b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 864ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 8658f99978dSLois Curfman McInnes 866b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 867b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 868b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 869ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 870ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 871ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 872ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 873b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8748f99978dSLois Curfman McInnes 875c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 876c8dd0c0dSLois Curfman McInnes 877c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 878c8dd0c0dSLois Curfman McInnes @*/ 8798f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 880c8dd0c0dSLois Curfman McInnes { 8818f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 882c8dd0c0dSLois Curfman McInnes 883c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 884c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 885c8dd0c0dSLois Curfman McInnes if (f) { 886c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 887c8dd0c0dSLois Curfman McInnes } 888c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 889c8dd0c0dSLois Curfman McInnes } 890c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 891c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 892c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 893b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8948f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 895c8dd0c0dSLois Curfman McInnes { 896c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8976831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8986831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 899c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 900c8dd0c0dSLois Curfman McInnes } 901c8dd0c0dSLois Curfman McInnes EXTERN_C_END 902c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 90304d965bbSLois Curfman McInnes /* 9046831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 90582bf6240SBarry Smith 90604d965bbSLois Curfman McInnes Input Parameter: 90704d965bbSLois Curfman McInnes . snes - the SNES context 90804d965bbSLois Curfman McInnes 90904d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 91004d965bbSLois Curfman McInnes */ 9115615d1e5SSatish Balay #undef __FUNC__ 912b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 913f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 9145e42470aSBarry Smith { 9156831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 916d132466eSBarry Smith int ierr; 9176daaf66cSBarry Smith 9183a40ed3dSBarry Smith PetscFunctionBegin; 9196831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 920d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 921d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 922d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 923d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 9243a40ed3dSBarry Smith PetscFunctionReturn(0); 9255e42470aSBarry Smith } 92604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92704d965bbSLois Curfman McInnes /* 9286831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 92904d965bbSLois Curfman McInnes 93004d965bbSLois Curfman McInnes Input Parameters: 93104d965bbSLois Curfman McInnes . SNES - the SNES context 93204d965bbSLois Curfman McInnes . viewer - visualization context 93304d965bbSLois Curfman McInnes 93404d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 93504d965bbSLois Curfman McInnes */ 9365615d1e5SSatish Balay #undef __FUNC__ 937b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 938e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 939a935fc98SLois Curfman McInnes { 9406831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 94119bcc07fSBarry Smith char *cstr; 94251695f54SLois Curfman McInnes int ierr; 9436831982aSBarry Smith PetscTruth isascii; 944a935fc98SLois Curfman McInnes 9453a40ed3dSBarry Smith PetscFunctionBegin; 9466831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 9470f5bd95cSBarry Smith if (isascii) { 94819bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 94919bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 95019bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 95119bcc07fSBarry Smith else cstr = "unknown"; 952d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 953d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9545cd90555SBarry Smith } else { 9550f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 95619bcc07fSBarry Smith } 9573a40ed3dSBarry Smith PetscFunctionReturn(0); 958a935fc98SLois Curfman McInnes } 95904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 96004d965bbSLois Curfman McInnes /* 9616831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 96204d965bbSLois Curfman McInnes 96304d965bbSLois Curfman McInnes Input Parameter: 96404d965bbSLois Curfman McInnes . snes - the SNES context 96504d965bbSLois Curfman McInnes 96604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 96704d965bbSLois Curfman McInnes */ 9685615d1e5SSatish Balay #undef __FUNC__ 969b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 970f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 9715e42470aSBarry Smith { 9726831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 9735e42470aSBarry Smith char ver[16]; 9745e42470aSBarry Smith double tmp; 975f1af5d2fSBarry Smith int ierr; 976f1af5d2fSBarry Smith PetscTruth flg; 9775e42470aSBarry Smith 9783a40ed3dSBarry Smith PetscFunctionBegin; 97909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 98025cdf11fSBarry Smith if (flg) { 9815e42470aSBarry Smith ls->alpha = tmp; 9825e42470aSBarry Smith } 98309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 98425cdf11fSBarry Smith if (flg) { 9855e42470aSBarry Smith ls->maxstep = tmp; 9865e42470aSBarry Smith } 98709d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 98825cdf11fSBarry Smith if (flg) { 9895e42470aSBarry Smith ls->steptol = tmp; 9905e42470aSBarry Smith } 99109d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 99225cdf11fSBarry Smith if (flg) { 993c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9940f5bd95cSBarry Smith 995c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 996c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 997c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 998c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9990f5bd95cSBarry Smith 10000f5bd95cSBarry Smith if (isbasic) { 1001af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 10020f5bd95cSBarry Smith } else if (isnonorms) { 1003af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 10040f5bd95cSBarry Smith } else if (isquad) { 1005af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 10060f5bd95cSBarry Smith } else if (iscubic) { 1007af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 10085e42470aSBarry Smith } 1009a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 10105e42470aSBarry Smith } 10113a40ed3dSBarry Smith PetscFunctionReturn(0); 10125e42470aSBarry Smith } 101304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 101404d965bbSLois Curfman McInnes /* 10156831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 10166831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 101704d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 101804d965bbSLois Curfman McInnes 101904d965bbSLois Curfman McInnes 102004d965bbSLois Curfman McInnes Input Parameter: 102104d965bbSLois Curfman McInnes . snes - the SNES context 102204d965bbSLois Curfman McInnes 102304d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 102404d965bbSLois Curfman McInnes */ 1025fb2e594dSBarry Smith EXTERN_C_BEGIN 10265615d1e5SSatish Balay #undef __FUNC__ 1027b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 1028f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 10295e42470aSBarry Smith { 103082bf6240SBarry Smith int ierr; 10316831982aSBarry Smith SNES_EQ_LS *neP; 10325e42470aSBarry Smith 10333a40ed3dSBarry Smith PetscFunctionBegin; 1034a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1035a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1036a8c6a408SBarry Smith } 103782bf6240SBarry Smith 1038f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 1039f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 1040f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 1041f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 1042f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 1043f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 1044f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 10455baf8537SBarry Smith snes->nwork = 0; 10465e42470aSBarry Smith 10476831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 10486831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 10495e42470aSBarry Smith snes->data = (void*)neP; 10505e42470aSBarry Smith neP->alpha = 1.e-4; 10515e42470aSBarry Smith neP->maxstep = 1.e8; 10525e42470aSBarry Smith neP->steptol = 1.e-12; 10535e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1054c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1055c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1056c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 105782bf6240SBarry Smith 10585d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10595d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 106082bf6240SBarry Smith 10613a40ed3dSBarry Smith PetscFunctionReturn(0); 10625e42470aSBarry Smith } 1063fb2e594dSBarry Smith EXTERN_C_END 106482bf6240SBarry Smith 106582bf6240SBarry Smith 106682bf6240SBarry Smith 1067