15e76c431SBarry Smith 270f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 35e42470aSBarry Smith 48a5d9ee4SBarry Smith /* 58a5d9ee4SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the function, 68a5d9ee4SBarry Smith but not a zero. In the case when one cannot compute J^T F we use the fact that 736669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 836669109SBarry Smith for this trick. 98a5d9ee4SBarry Smith */ 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 12dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 138a5d9ee4SBarry Smith { 148a5d9ee4SBarry Smith PetscReal a1; 15dfbe8321SBarry Smith PetscErrorCode ierr; 1636669109SBarry Smith PetscTruth hastranspose; 178a5d9ee4SBarry Smith 188a5d9ee4SBarry Smith PetscFunctionBegin; 198a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2036669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2136669109SBarry Smith if (hastranspose) { 228a5d9ee4SBarry Smith /* Compute || J^T F|| */ 238a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 248a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 254b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm); 268a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2736669109SBarry Smith } else { 2836669109SBarry Smith Vec work; 29ea709b57SSatish Balay PetscScalar result; 3036669109SBarry Smith PetscReal wnorm; 3136669109SBarry Smith 3236669109SBarry Smith ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr); 3336669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3536669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 3736669109SBarry Smith ierr = VecDestroy(work);CHKERRQ(ierr); 3836669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 394b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1); 4036669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4136669109SBarry Smith } 428a5d9ee4SBarry Smith PetscFunctionReturn(0); 438a5d9ee4SBarry Smith } 448a5d9ee4SBarry Smith 4574637425SBarry Smith /* 465ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4774637425SBarry Smith */ 484a2ae208SSatish Balay #undef __FUNCT__ 494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 50dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5174637425SBarry Smith { 5274637425SBarry Smith PetscReal a1,a2; 53dfbe8321SBarry Smith PetscErrorCode ierr; 5474637425SBarry Smith PetscTruth hastranspose; 55ea709b57SSatish Balay PetscScalar mone = -1.0; 5674637425SBarry Smith 5774637425SBarry Smith PetscFunctionBegin; 5874637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 5974637425SBarry Smith if (hastranspose) { 6074637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6174637425SBarry Smith ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr); 6274637425SBarry Smith 6374637425SBarry Smith /* Compute || J^T W|| */ 6474637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6574637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 679994e62eSSatish Balay if (a1 != 0) { 684b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1); 6985471664SBarry Smith } 7074637425SBarry Smith } 7174637425SBarry Smith PetscFunctionReturn(0); 7274637425SBarry Smith } 7374637425SBarry Smith 7404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7504d965bbSLois Curfman McInnes 7604d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7794b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7804d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7904d965bbSLois Curfman McInnes respectively. 8004d965bbSLois Curfman McInnes 81fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8204d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8304d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 84fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8504d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8604d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 874b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 884b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 8904d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9004d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9104d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9204d965bbSLois Curfman McInnes for all nonlinear solvers. 9304d965bbSLois Curfman McInnes 9404d965bbSLois Curfman McInnes Another key routine is: 9504d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9604d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9704d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9804d965bbSLois Curfman McInnes SNESSolve() if necessary. 9904d965bbSLois Curfman McInnes 10004d965bbSLois Curfman McInnes Additional basic routines are: 10104d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10204d965bbSLois Curfman McInnes have actually been used. 10304d965bbSLois Curfman McInnes These are called by application codes via the interface routines 104186905e3SBarry Smith SNESView(). 10504d965bbSLois Curfman McInnes 10604d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10704d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10804d965bbSLois Curfman McInnes above description applies to these categories also. 10904d965bbSLois Curfman McInnes 11004d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1115e42470aSBarry Smith /* 1124b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11304d965bbSLois Curfman McInnes method with a line search. 1145e76c431SBarry Smith 11504d965bbSLois Curfman McInnes Input Parameters: 11604d965bbSLois Curfman McInnes . snes - the SNES context 1175e76c431SBarry Smith 11804d965bbSLois Curfman McInnes Output Parameter: 119c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12004d965bbSLois Curfman McInnes 12104d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1225e76c431SBarry Smith 1235e76c431SBarry Smith Notes: 1245e76c431SBarry Smith This implements essentially a truncated Newton method with a 1255e76c431SBarry Smith line search. By default a cubic backtracking line search 1265e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1275e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 128393d2d9aSLois Curfman McInnes and Schnabel. 1295e76c431SBarry Smith */ 1304a2ae208SSatish Balay #undef __FUNCT__ 1314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 132dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1335e76c431SBarry Smith { 1344b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1356849ba73SBarry Smith PetscErrorCode ierr; 136a7cc72afSBarry Smith PetscInt maxits,i,lits; 137a7cc72afSBarry Smith PetscTruth lssucceed; 138112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1405e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 141c293cc10SBarry Smith KSP ksp; 1425e76c431SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 14494b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 145184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 146184914b5SBarry Smith 1475e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1485e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1505e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1515e42470aSBarry Smith G = snes->work[1]; 1525e42470aSBarry Smith W = snes->work[2]; 1535e76c431SBarry Smith 1544c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1554c49b128SBarry Smith snes->iter = 0; 1564c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1575334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 158cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1590f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1605e42470aSBarry Smith snes->norm = fnorm; 1610f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16284c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16394a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1645e76c431SBarry Smith 16570441072SBarry Smith if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16694a9d846SBarry Smith 167d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 168d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 169d034289bSBarry Smith 1705e76c431SBarry Smith for (i=0; i<maxits; i++) { 1715e76c431SBarry Smith 172329e7e40SMatthew Knepley /* Call general purpose update function */ 173abc0a331SBarry Smith if (snes->update) { 174329e7e40SMatthew Knepley ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr); 175de8cb200SMatthew Knepley } 176329e7e40SMatthew Knepley 177ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1785334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 17994b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 18023ce1328SBarry Smith ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 181c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 18274637425SBarry Smith 183b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 18474637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 18585471664SBarry Smith } 18674637425SBarry Smith 18790cbd584SBarry Smith /* should check what happened to the linear solve? */ 1883505fcc1SBarry Smith snes->linear_its += lits; 18977431f27SBarry Smith PetscLogInfo(snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits); 190ea4d3ed3SLois Curfman McInnes 191ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 192ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 193ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 194ea4d3ed3SLois Curfman McInnes */ 19581b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 196a7cc72afSBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 197a7cc72afSBarry Smith PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed); 198a3b891d8SBarry Smith 199a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 200a3b891d8SBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 201a3b891d8SBarry Smith fnorm = gnorm; 202a3b891d8SBarry Smith 2035ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2045ed2d596SBarry Smith snes->iter = i+1; 2055ed2d596SBarry Smith snes->norm = fnorm; 2065ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2075ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2085ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2095ed2d596SBarry Smith 210a7cc72afSBarry Smith if (!lssucceed) { 2118a5d9ee4SBarry Smith PetscTruth ismin; 21250ffb88aSMatthew Knepley 21350ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2143505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2158a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2168a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2173505fcc1SBarry Smith break; 2183505fcc1SBarry Smith } 21950ffb88aSMatthew Knepley } 2205e76c431SBarry Smith 2215e76c431SBarry Smith /* Test for convergence */ 22229e0b56fSBarry Smith if (snes->converged) { 22329e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2243505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2253505fcc1SBarry Smith if (snes->reason) { 22616c95cacSBarry Smith break; 22716c95cacSBarry Smith } 22816c95cacSBarry Smith } 22929e0b56fSBarry Smith } 23039e2f89bSBarry Smith if (X != snes->vec_sol) { 231393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 232e15f7bb5SBarry Smith } 233e15f7bb5SBarry Smith if (F != snes->vec_func) { 234e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 235e15f7bb5SBarry Smith } 23639e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 23739e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 23852392280SLois Curfman McInnes if (i == maxits) { 23977431f27SBarry Smith PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits); 2403505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 24152392280SLois Curfman McInnes } 2423a40ed3dSBarry Smith PetscFunctionReturn(0); 2435e76c431SBarry Smith } 24404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 24504d965bbSLois Curfman McInnes /* 2464b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2474b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 24804d965bbSLois Curfman McInnes 24904d965bbSLois Curfman McInnes Input Parameter: 25004d965bbSLois Curfman McInnes . snes - the SNES context 25104d965bbSLois Curfman McInnes . x - the solution vector 25204d965bbSLois Curfman McInnes 25304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 25404d965bbSLois Curfman McInnes 25504d965bbSLois Curfman McInnes Notes: 2564b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 25704d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 25804d965bbSLois Curfman McInnes the call to SNESSolve(). 25904d965bbSLois Curfman McInnes */ 2604a2ae208SSatish Balay #undef __FUNCT__ 2614b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 262dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2635e76c431SBarry Smith { 264dfbe8321SBarry Smith PetscErrorCode ierr; 2653a40ed3dSBarry Smith 2663a40ed3dSBarry Smith PetscFunctionBegin; 26781b6cf68SLois Curfman McInnes snes->nwork = 4; 268d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 269b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 27081b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2725e76c431SBarry Smith } 27304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27404d965bbSLois Curfman McInnes /* 2754b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2764b27c08aSLois Curfman McInnes with SNESCreate_LS(). 27704d965bbSLois Curfman McInnes 27804d965bbSLois Curfman McInnes Input Parameter: 27904d965bbSLois Curfman McInnes . snes - the SNES context 28004d965bbSLois Curfman McInnes 281de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 28204d965bbSLois Curfman McInnes */ 2834a2ae208SSatish Balay #undef __FUNCT__ 2844b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 285dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 2865e76c431SBarry Smith { 287dfbe8321SBarry Smith PetscErrorCode ierr; 2883a40ed3dSBarry Smith 2893a40ed3dSBarry Smith PetscFunctionBegin; 2905baf8537SBarry Smith if (snes->nwork) { 2914b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 29221c89e3eSBarry Smith } 2935c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 2955e76c431SBarry Smith } 29604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2974a2ae208SSatish Balay #undef __FUNCT__ 2984a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch" 29982bf6240SBarry Smith 3004b828684SBarry Smith /*@C 3015e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 3025e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3035e76c431SBarry Smith to serve as a template and is not recommended for general use. 3045e76c431SBarry Smith 305c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 306c7afd0dbSLois Curfman McInnes 3075e76c431SBarry Smith Input Parameters: 308c7afd0dbSLois Curfman McInnes + snes - nonlinear context 309af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3105e76c431SBarry Smith . x - current iterate 3115e76c431SBarry Smith . f - residual evaluated at x 3125e76c431SBarry Smith . y - search direction (contains new iterate on output) 3135e76c431SBarry Smith . w - work vector 314c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3155e76c431SBarry Smith 316c4a48953SLois Curfman McInnes Output Parameters: 317c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3185e76c431SBarry Smith . y - new iterate (contains search direction on input) 3195e76c431SBarry Smith . gnorm - 2-norm of g 3205e76c431SBarry Smith . ynorm - 2-norm of search length 321a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 322fee21e36SBarry Smith 323c4a48953SLois Curfman McInnes Options Database Key: 3244b27c08aSLois Curfman McInnes . -snes_ls basic - Activates SNESNoLineSearch() 325c4a48953SLois Curfman McInnes 32636851e7fSLois Curfman McInnes Level: advanced 32736851e7fSLois Curfman McInnes 32828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 32928ae5a21SLois Curfman McInnes 330f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 331af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3325e76c431SBarry Smith @*/ 333a7cc72afSBarry Smith PetscErrorCode SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 3345e76c431SBarry Smith { 335dfbe8321SBarry Smith PetscErrorCode ierr; 336ea709b57SSatish Balay PetscScalar mone = -1.0; 3374b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3388f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3395334005bSBarry Smith 3403a40ed3dSBarry Smith PetscFunctionBegin; 341a7cc72afSBarry Smith *flag = PETSC_TRUE; 342d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 343a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 344ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3458f99978dSLois Curfman McInnes if (neP->CheckStep) { 3468f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3478f99978dSLois Curfman McInnes } 348ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 349a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 350d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3513a40ed3dSBarry Smith PetscFunctionReturn(0); 3525e76c431SBarry Smith } 35304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3544a2ae208SSatish Balay #undef __FUNCT__ 3554a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms" 35682bf6240SBarry Smith 35729e0b56fSBarry Smith /*@C 35829e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 35929e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 36029e0b56fSBarry Smith even compute the norm of the function or search direction; this 36129e0b56fSBarry Smith is intended only when you know the full step is fine and are 36229e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 363c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 364c7afd0dbSLois Curfman McInnes 365c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 36629e0b56fSBarry Smith 36729e0b56fSBarry Smith Input Parameters: 368c7afd0dbSLois Curfman McInnes + snes - nonlinear context 369af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 37029e0b56fSBarry Smith . x - current iterate 37129e0b56fSBarry Smith . f - residual evaluated at x 37229e0b56fSBarry Smith . y - search direction (contains new iterate on output) 37329e0b56fSBarry Smith . w - work vector 374c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 37529e0b56fSBarry Smith 37629e0b56fSBarry Smith Output Parameters: 377c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 37829e0b56fSBarry Smith . gnorm - not changed 37929e0b56fSBarry Smith . ynorm - not changed 380a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 381fee21e36SBarry Smith 38229e0b56fSBarry Smith Options Database Key: 3834b27c08aSLois Curfman McInnes . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 38429e0b56fSBarry Smith 3858cbba510SBarry Smith Notes: 386ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 387ea56f5baSLois Curfman McInnes either the options 388ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 389ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 390329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 391329f5518SBarry Smith otherwise, the SNES solver will generate an error. 392329f5518SBarry Smith 393329f5518SBarry Smith During the final iteration this will not evaluate the function at 394329f5518SBarry Smith the solution point. This is to save a function evaluation while 395329f5518SBarry Smith using pseudo-timestepping. 3968cbba510SBarry Smith 397ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 398ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 399ea56f5baSLois Curfman McInnes correct, since they are not computed. 400ea56f5baSLois Curfman McInnes 401ea56f5baSLois Curfman McInnes Level: advanced 4028cbba510SBarry Smith 40329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 40429e0b56fSBarry Smith 40529e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 40629e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 40729e0b56fSBarry Smith @*/ 408a7cc72afSBarry Smith PetscErrorCode SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 40929e0b56fSBarry Smith { 410dfbe8321SBarry Smith PetscErrorCode ierr; 411ea709b57SSatish Balay PetscScalar mone = -1.0; 4124b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4138f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 41429e0b56fSBarry Smith 4153a40ed3dSBarry Smith PetscFunctionBegin; 416a7cc72afSBarry Smith *flag = PETSC_TRUE; 417d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 41829e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4198f99978dSLois Curfman McInnes if (neP->CheckStep) { 4208f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4218f99978dSLois Curfman McInnes } 422329f5518SBarry Smith 423329f5518SBarry Smith /* don't evaluate function the last time through */ 424329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 42529e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 426329f5518SBarry Smith } 427d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4283a40ed3dSBarry Smith PetscFunctionReturn(0); 42929e0b56fSBarry Smith } 43004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4314a2ae208SSatish Balay #undef __FUNCT__ 4324a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch" 4334b828684SBarry Smith /*@C 434f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4355e76c431SBarry Smith 436c7afd0dbSLois Curfman McInnes Collective on SNES 437c7afd0dbSLois Curfman McInnes 4385e76c431SBarry Smith Input Parameters: 439c7afd0dbSLois Curfman McInnes + snes - nonlinear context 440af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4415e76c431SBarry Smith . x - current iterate 4425e76c431SBarry Smith . f - residual evaluated at x 4435e76c431SBarry Smith . y - search direction (contains new iterate on output) 4445e76c431SBarry Smith . w - work vector 445c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4465e76c431SBarry Smith 447393d2d9aSLois Curfman McInnes Output Parameters: 448c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4495e76c431SBarry Smith . y - new iterate (contains search direction on input) 4505e76c431SBarry Smith . gnorm - 2-norm of g 4515e76c431SBarry Smith . ynorm - 2-norm of search length 452a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 453fee21e36SBarry Smith 454c4a48953SLois Curfman McInnes Options Database Key: 4554b27c08aSLois Curfman McInnes $ -snes_ls cubic - Activates SNESCubicLineSearch() 456c4a48953SLois Curfman McInnes 4575e76c431SBarry Smith Notes: 4585e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4595e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4605e76c431SBarry Smith 46136851e7fSLois Curfman McInnes Level: advanced 46236851e7fSLois Curfman McInnes 46328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 46428ae5a21SLois Curfman McInnes 465af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4665e76c431SBarry Smith @*/ 467a7cc72afSBarry Smith PetscErrorCode SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 4685e76c431SBarry Smith { 469406556e6SLois Curfman McInnes /* 470406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 471406556e6SLois Curfman McInnes minimization problem: 472406556e6SLois Curfman McInnes min z(x): R^n -> R, 473406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 474406556e6SLois Curfman McInnes */ 475406556e6SLois Curfman McInnes 4765ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 477329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 478aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 479ea709b57SSatish Balay PetscScalar cinitslope,clambda; 4806b5873e3SBarry Smith #endif 481dfbe8321SBarry Smith PetscErrorCode ierr; 482a7cc72afSBarry Smith PetscInt count; 4834b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 484ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 4858f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4865e76c431SBarry Smith 4873a40ed3dSBarry Smith PetscFunctionBegin; 488d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 489a7cc72afSBarry Smith *flag = PETSC_TRUE; 4905e76c431SBarry Smith alpha = neP->alpha; 4915e76c431SBarry Smith maxstep = neP->maxstep; 4925e76c431SBarry Smith steptol = neP->steptol; 4935e76c431SBarry Smith 494cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 49563b9fa5eSBarry Smith if (*ynorm == 0.0) { 49663b9fa5eSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n"); 497a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 498a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 499a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 500a7cc72afSBarry Smith *flag = PETSC_FALSE; 501ad922ac9SBarry Smith goto theend1; 50294a9d846SBarry Smith } 5035e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5045e42470aSBarry Smith scale = maxstep/(*ynorm); 505aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 506b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 5076b5873e3SBarry Smith #else 508b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 5096b5873e3SBarry Smith #endif 510393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5115e76c431SBarry Smith *ynorm = maxstep; 5125e76c431SBarry Smith } 513ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5145ca10a37SBarry Smith minlambda = steptol/rellength; 515a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 516aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 517a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 518329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5195e42470aSBarry Smith #else 520a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5215e42470aSBarry Smith #endif 5225e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5235e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5245e76c431SBarry Smith 525393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5265334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 52778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 528313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5295d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 530393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 531b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 532ad922ac9SBarry Smith goto theend1; 5335e76c431SBarry Smith } 5345e76c431SBarry Smith 5355e76c431SBarry Smith /* Fit points with quadratic */ 536313b5538SBarry Smith lambda = 1.0; 537a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5385e76c431SBarry Smith lambdaprev = lambda; 5395e76c431SBarry Smith gnormprev = *gnorm; 540329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 541ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 542ddd12767SBarry Smith else lambda = lambdatemp; 543393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 544ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 545aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 546ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5475e42470aSBarry Smith #else 548ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5495e42470aSBarry Smith #endif 55078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 551cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5525ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 553393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 5545ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda); 555ad922ac9SBarry Smith goto theend1; 5565e76c431SBarry Smith } 5575e76c431SBarry Smith 5585e76c431SBarry Smith /* Fit points with cubic */ 5595e76c431SBarry Smith count = 1; 5608229bfc2SMatthew Knepley while (count < 10000) { 5615e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 56277431f27SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count); 5635ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope); 564f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 565a7cc72afSBarry Smith *flag = PETSC_FALSE; break; 5665e76c431SBarry Smith } 5675d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 5685d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 569ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5702b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5715e76c431SBarry Smith d = b*b - 3*a*initslope; 5725e76c431SBarry Smith if (d < 0.0) d = 0.0; 5735e76c431SBarry Smith if (a == 0.0) { 5745e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5755e76c431SBarry Smith } else { 5765e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5775e76c431SBarry Smith } 5785e76c431SBarry Smith lambdaprev = lambda; 5795e76c431SBarry Smith gnormprev = *gnorm; 580329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 581329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5825e76c431SBarry Smith else lambda = lambdatemp; 583393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 584ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 585aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 586ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 587393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5885e42470aSBarry Smith #else 589ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5905e42470aSBarry Smith #endif 59178b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 592cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5935ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 594393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 5955ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda); 5962715565aSLois Curfman McInnes goto theend1; 5972b022350SLois Curfman McInnes } else { 5985ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda); 5995e76c431SBarry Smith } 6005e76c431SBarry Smith count++; 6015e76c431SBarry Smith } 6028229bfc2SMatthew Knepley if (count >= 10000) { 603cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6048229bfc2SMatthew Knepley } 605ad922ac9SBarry Smith theend1: 6068f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6078f99978dSLois Curfman McInnes if (neP->CheckStep) { 6088f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 609abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 6108f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6118f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6128f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6138f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6148f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6158f99978dSLois Curfman McInnes } 6168f99978dSLois Curfman McInnes } 617d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6183a40ed3dSBarry Smith PetscFunctionReturn(0); 6195e76c431SBarry Smith } 62004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6214a2ae208SSatish Balay #undef __FUNCT__ 6224a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch" 6234b828684SBarry Smith /*@C 624f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6255e76c431SBarry Smith 626c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 627c7afd0dbSLois Curfman McInnes 6285e42470aSBarry Smith Input Parameters: 629c7afd0dbSLois Curfman McInnes + snes - the SNES context 630af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6315e42470aSBarry Smith . x - current iterate 6325e42470aSBarry Smith . f - residual evaluated at x 6335e42470aSBarry Smith . y - search direction (contains new iterate on output) 6345e42470aSBarry Smith . w - work vector 635c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6365e42470aSBarry Smith 637c4a48953SLois Curfman McInnes Output Parameters: 638c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6395e42470aSBarry Smith . y - new iterate (contains search direction on input) 6405e42470aSBarry Smith . gnorm - 2-norm of g 6415e42470aSBarry Smith . ynorm - 2-norm of search length 642a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 643fee21e36SBarry Smith 644c4a48953SLois Curfman McInnes Options Database Key: 6454b27c08aSLois Curfman McInnes . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 646c4a48953SLois Curfman McInnes 6475e42470aSBarry Smith Notes: 6484b27c08aSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNESLS method. 6495e42470aSBarry Smith 65036851e7fSLois Curfman McInnes Level: advanced 65136851e7fSLois Curfman McInnes 652f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 653f59ffdedSLois Curfman McInnes 654af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6555e42470aSBarry Smith @*/ 656a7cc72afSBarry Smith PetscErrorCode SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 6575e76c431SBarry Smith { 658406556e6SLois Curfman McInnes /* 659406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 660406556e6SLois Curfman McInnes minimization problem: 661406556e6SLois Curfman McInnes min z(x): R^n -> R, 662406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 663406556e6SLois Curfman McInnes */ 6645ca10a37SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 665aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 666ea709b57SSatish Balay PetscScalar cinitslope,clambda; 6676b5873e3SBarry Smith #endif 668dfbe8321SBarry Smith PetscErrorCode ierr; 669a7cc72afSBarry Smith PetscInt count; 6704b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 671ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 6728f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6735e76c431SBarry Smith 6743a40ed3dSBarry Smith PetscFunctionBegin; 675d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 676a7cc72afSBarry Smith *flag = PETSC_TRUE; 6775e42470aSBarry Smith alpha = neP->alpha; 6785e42470aSBarry Smith maxstep = neP->maxstep; 6795e42470aSBarry Smith steptol = neP->steptol; 6805e76c431SBarry Smith 6813505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 68263b9fa5eSBarry Smith if (*ynorm == 0.0) { 683b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 684b37302c6SLois Curfman McInnes *gnorm = fnorm; 685b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 686b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 687a7cc72afSBarry Smith *flag = PETSC_FALSE; 688ad922ac9SBarry Smith goto theend2; 68994a9d846SBarry Smith } 6905e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6915e42470aSBarry Smith scale = maxstep/(*ynorm); 692393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6935e42470aSBarry Smith *ynorm = maxstep; 6945e76c431SBarry Smith } 695ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 6965ca10a37SBarry Smith minlambda = steptol/rellength; 697a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 698aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 699a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 700329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7015e42470aSBarry Smith #else 702a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7035e42470aSBarry Smith #endif 7045e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7055e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7065e42470aSBarry Smith 707393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 7085334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 70978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 710cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7115d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 712393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 713b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 714ad922ac9SBarry Smith goto theend2; 7155e42470aSBarry Smith } 7165e42470aSBarry Smith 7175e42470aSBarry Smith /* Fit points with quadratic */ 718313b5538SBarry Smith lambda = 1.0; 7195e42470aSBarry Smith count = 1; 7205ca10a37SBarry Smith while (PETSC_TRUE) { 7215e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 72277431f27SBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count); 723b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 724f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 725a7cc72afSBarry Smith *flag = PETSC_FALSE; break; 7265e42470aSBarry Smith } 727a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 728329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 729329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 730329f5518SBarry Smith else lambda = lambdatemp; 731393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7323505fcc1SBarry Smith lambdaneg = -lambda; 733aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7343505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7355e42470aSBarry Smith #else 7363505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7375e42470aSBarry Smith #endif 73878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 739cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7405ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 741393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 742b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 74306259719SBarry Smith break; 7445e42470aSBarry Smith } 7455e42470aSBarry Smith count++; 7465e42470aSBarry Smith } 747ad922ac9SBarry Smith theend2: 7488f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7498f99978dSLois Curfman McInnes if (neP->CheckStep) { 7508f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 751abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 7528f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7538f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7548f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7558f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7568f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7578f99978dSLois Curfman McInnes } 7588f99978dSLois Curfman McInnes } 759d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7603a40ed3dSBarry Smith PetscFunctionReturn(0); 7615e76c431SBarry Smith } 7622343ba6eSBarry Smith 76304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7644a2ae208SSatish Balay #undef __FUNCT__ 7654a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 766c9e6a524SLois Curfman McInnes /*@C 76777c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7684b27c08aSLois Curfman McInnes by the method SNESLS. 7695e76c431SBarry Smith 7705e76c431SBarry Smith Input Parameters: 771c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 772af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 773c7afd0dbSLois Curfman McInnes - func - pointer to int function 7745e76c431SBarry Smith 775fee21e36SBarry Smith Collective on SNES 776fee21e36SBarry Smith 777c4a48953SLois Curfman McInnes Available Routines: 778c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 779f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 780f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 781af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7825e76c431SBarry Smith 783c4a48953SLois Curfman McInnes Options Database Keys: 7844b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 7854b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 7864b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 7874b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 7883304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 7893304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 790c4a48953SLois Curfman McInnes 7915e76c431SBarry Smith Calling sequence of func: 792bd208895SLois Curfman McInnes .vb 793af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 794a7cc72afSBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 795bd208895SLois Curfman McInnes .ve 7965e76c431SBarry Smith 7975e76c431SBarry Smith Input parameters for func: 798c7afd0dbSLois Curfman McInnes + snes - nonlinear context 799af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8005e76c431SBarry Smith . x - current iterate 8015e76c431SBarry Smith . f - residual evaluated at x 8025e76c431SBarry Smith . y - search direction (contains new iterate on output) 8035e76c431SBarry Smith . w - work vector 804c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 8055e76c431SBarry Smith 8065e76c431SBarry Smith Output parameters for func: 807c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 8085e76c431SBarry Smith . y - new iterate (contains search direction on input) 8095e76c431SBarry Smith . gnorm - 2-norm of g 8105e76c431SBarry Smith . ynorm - 2-norm of search length 811a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 812f59ffdedSLois Curfman McInnes 81336851e7fSLois Curfman McInnes Level: advanced 81436851e7fSLois Curfman McInnes 815f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 816f59ffdedSLois Curfman McInnes 8175d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 8185d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 8195e76c431SBarry Smith @*/ 820a7cc72afSBarry Smith PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 8215e76c431SBarry Smith { 822a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 82382bf6240SBarry Smith 8243a40ed3dSBarry Smith PetscFunctionBegin; 825c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 82682bf6240SBarry Smith if (f) { 827af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 82882bf6240SBarry Smith } 8293a40ed3dSBarry Smith PetscFunctionReturn(0); 8305e76c431SBarry Smith } 8318e019c35SBarry Smith 832a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 83304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 834fb2e594dSBarry Smith EXTERN_C_BEGIN 8354a2ae208SSatish Balay #undef __FUNCT__ 8364a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS" 837dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx) 83882bf6240SBarry Smith { 83982bf6240SBarry Smith PetscFunctionBegin; 8404b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 8414b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 84282bf6240SBarry Smith PetscFunctionReturn(0); 84382bf6240SBarry Smith } 844fb2e594dSBarry Smith EXTERN_C_END 84504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8464a2ae208SSatish Balay #undef __FUNCT__ 8474a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck" 848c8dd0c0dSLois Curfman McInnes /*@C 849530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 8504b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 851c8dd0c0dSLois Curfman McInnes 852c8dd0c0dSLois Curfman McInnes Input Parameters: 853c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 854c8dd0c0dSLois Curfman McInnes . func - pointer to int function 855c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 856c8dd0c0dSLois Curfman McInnes 857c8dd0c0dSLois Curfman McInnes Collective on SNES 858c8dd0c0dSLois Curfman McInnes 859c8dd0c0dSLois Curfman McInnes Calling sequence of func: 860c8dd0c0dSLois Curfman McInnes .vb 861b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 862c8dd0c0dSLois Curfman McInnes .ve 863b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 864b1ae27eaSLois Curfman McInnes on failure. 865c8dd0c0dSLois Curfman McInnes 866c8dd0c0dSLois Curfman McInnes Input parameters for func: 867c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 868c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8698f99978dSLois Curfman McInnes - x - current candidate iterate 870c8dd0c0dSLois Curfman McInnes 871c8dd0c0dSLois Curfman McInnes Output parameters for func: 872c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 873c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 874c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 875c8dd0c0dSLois Curfman McInnes 876c8dd0c0dSLois Curfman McInnes Level: advanced 877c8dd0c0dSLois Curfman McInnes 8788f99978dSLois Curfman McInnes Notes: 879b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 880b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 881b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 882b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 883ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 8848f99978dSLois Curfman McInnes 885b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 886b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 887b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 888ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 889ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 890ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 891ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 892b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8938f99978dSLois Curfman McInnes 894c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 895c8dd0c0dSLois Curfman McInnes 896c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 897c8dd0c0dSLois Curfman McInnes @*/ 8986849ba73SBarry Smith PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 899c8dd0c0dSLois Curfman McInnes { 9006849ba73SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*); 901c8dd0c0dSLois Curfman McInnes 902c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 903c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 904c8dd0c0dSLois Curfman McInnes if (f) { 905c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 906c8dd0c0dSLois Curfman McInnes } 907c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 908c8dd0c0dSLois Curfman McInnes } 909c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9106849ba73SBarry Smith typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/ 911c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 9124a2ae208SSatish Balay #undef __FUNCT__ 9134a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS" 914dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx) 915c8dd0c0dSLois Curfman McInnes { 916c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 9174b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 9184b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 919c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 920c8dd0c0dSLois Curfman McInnes } 921c8dd0c0dSLois Curfman McInnes EXTERN_C_END 922c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92304d965bbSLois Curfman McInnes /* 9244b27c08aSLois Curfman McInnes SNESPrintHelp_LS - Prints all options for the SNES_LS method. 925329e7e40SMatthew Knepley 926329e7e40SMatthew Knepley Input Parameter: 927329e7e40SMatthew Knepley . snes - the SNES context 928329e7e40SMatthew Knepley 929329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 930329e7e40SMatthew Knepley */ 931329e7e40SMatthew Knepley #undef __FUNCT__ 9324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS" 9336849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p) 934329e7e40SMatthew Knepley { 9354b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 936329e7e40SMatthew Knepley 937329e7e40SMatthew Knepley PetscFunctionBegin; 9384b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 9394b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 9404b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 9414b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 9424b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 943329e7e40SMatthew Knepley PetscFunctionReturn(0); 944329e7e40SMatthew Knepley } 945329e7e40SMatthew Knepley 946329e7e40SMatthew Knepley /* 9474b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 94804d965bbSLois Curfman McInnes 94904d965bbSLois Curfman McInnes Input Parameters: 95004d965bbSLois Curfman McInnes . SNES - the SNES context 95104d965bbSLois Curfman McInnes . viewer - visualization context 95204d965bbSLois Curfman McInnes 95304d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 95404d965bbSLois Curfman McInnes */ 9554a2ae208SSatish Balay #undef __FUNCT__ 9564b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 9576849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 958a935fc98SLois Curfman McInnes { 9594b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 9602fc52814SBarry Smith const char *cstr; 961dfbe8321SBarry Smith PetscErrorCode ierr; 96232077d6dSBarry Smith PetscTruth iascii; 963a935fc98SLois Curfman McInnes 9643a40ed3dSBarry Smith PetscFunctionBegin; 96532077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 96632077d6dSBarry Smith if (iascii) { 96719bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 96819bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 96919bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 97019bcc07fSBarry Smith else cstr = "unknown"; 971b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 972b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9735cd90555SBarry Smith } else { 97479a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 97519bcc07fSBarry Smith } 9763a40ed3dSBarry Smith PetscFunctionReturn(0); 977a935fc98SLois Curfman McInnes } 97804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 97904d965bbSLois Curfman McInnes /* 9804b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 98104d965bbSLois Curfman McInnes 98204d965bbSLois Curfman McInnes Input Parameter: 98304d965bbSLois Curfman McInnes . snes - the SNES context 98404d965bbSLois Curfman McInnes 98504d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 98604d965bbSLois Curfman McInnes */ 9874a2ae208SSatish Balay #undef __FUNCT__ 9884b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 9896849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 9905e42470aSBarry Smith { 9914b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 992e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 993dfbe8321SBarry Smith PetscErrorCode ierr; 994a7cc72afSBarry Smith PetscInt indx; 995f1af5d2fSBarry Smith PetscTruth flg; 9965e42470aSBarry Smith 9973a40ed3dSBarry Smith PetscFunctionBegin; 998b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 9994b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 10004b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 10014b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1002186905e3SBarry Smith 100322e36657SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 100425cdf11fSBarry Smith if (flg) { 100522e36657SBarry Smith switch (indx) { 1006b49fd9e1SBarry Smith case 0: 1007af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1008b49fd9e1SBarry Smith break; 1009b49fd9e1SBarry Smith case 1: 1010af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1011b49fd9e1SBarry Smith break; 1012b49fd9e1SBarry Smith case 2: 1013af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1014b49fd9e1SBarry Smith break; 1015b49fd9e1SBarry Smith case 3: 1016af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1017b49fd9e1SBarry Smith break; 10185e42470aSBarry Smith } 10195e42470aSBarry Smith } 1020b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10213a40ed3dSBarry Smith PetscFunctionReturn(0); 10225e42470aSBarry Smith } 102304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1024ebe3b25bSBarry Smith /*MC 1025ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 102604d965bbSLois Curfman McInnes 1027ebe3b25bSBarry Smith Options Database: 1028ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1029ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1030ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1031ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1032ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1033ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 103404d965bbSLois Curfman McInnes 1035ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 103604d965bbSLois Curfman McInnes 1037ee3001cbSBarry Smith Level: beginner 1038ee3001cbSBarry Smith 1039ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1040ebe3b25bSBarry Smith SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1041ebe3b25bSBarry Smith SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1042ebe3b25bSBarry Smith 1043ebe3b25bSBarry Smith M*/ 1044fb2e594dSBarry Smith EXTERN_C_BEGIN 10454a2ae208SSatish Balay #undef __FUNCT__ 10464b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 1047dfbe8321SBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 10485e42470aSBarry Smith { 1049dfbe8321SBarry Smith PetscErrorCode ierr; 10504b27c08aSLois Curfman McInnes SNES_LS *neP; 10515e42470aSBarry Smith 10523a40ed3dSBarry Smith PetscFunctionBegin; 10534b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 10544b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 10554b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 10564b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 10574b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 10584b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 10594b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 10605baf8537SBarry Smith snes->nwork = 0; 10615e42470aSBarry Smith 10624b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 1063*52e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 10645e42470aSBarry Smith snes->data = (void*)neP; 10655e42470aSBarry Smith neP->alpha = 1.e-4; 10665e42470aSBarry Smith neP->maxstep = 1.e8; 10675e42470aSBarry Smith neP->steptol = 1.e-12; 10685e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1069c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1070c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1071c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 107282bf6240SBarry Smith 10735d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10745d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 107582bf6240SBarry Smith 10763a40ed3dSBarry Smith PetscFunctionReturn(0); 10775e42470aSBarry Smith } 1078fb2e594dSBarry Smith EXTERN_C_END 107982bf6240SBarry Smith 108082bf6240SBarry Smith 108182bf6240SBarry Smith 1082