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); 269efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 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 591*ed98166eSMatthew Knepley if (snes->nfuncs > snes->max_funcs) { 592*ed98166eSMatthew Knepley PetscLogInfo(snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, but unable to find good step length! %D \n",count); 593*ed98166eSMatthew Knepley 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); 594*ed98166eSMatthew Knepley ierr = VecCopy(x,y);CHKERRQ(ierr); 595*ed98166eSMatthew Knepley *flag = PETSC_FALSE; 596*ed98166eSMatthew Knepley break; 597*ed98166eSMatthew Knepley } else { 59878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 599*ed98166eSMatthew Knepley } 600cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 6015ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 602393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 6035ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda); 6042715565aSLois Curfman McInnes goto theend1; 6052b022350SLois Curfman McInnes } else { 6065ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda); 6075e76c431SBarry Smith } 6085e76c431SBarry Smith count++; 6095e76c431SBarry Smith } 6108229bfc2SMatthew Knepley if (count >= 10000) { 611cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6128229bfc2SMatthew Knepley } 613ad922ac9SBarry Smith theend1: 6148f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6158f99978dSLois Curfman McInnes if (neP->CheckStep) { 6168f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 617abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 6188f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6198f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6208f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6218f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6228f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6238f99978dSLois Curfman McInnes } 6248f99978dSLois Curfman McInnes } 625d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6263a40ed3dSBarry Smith PetscFunctionReturn(0); 6275e76c431SBarry Smith } 62804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6294a2ae208SSatish Balay #undef __FUNCT__ 6304a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch" 6314b828684SBarry Smith /*@C 632f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6335e76c431SBarry Smith 634c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 635c7afd0dbSLois Curfman McInnes 6365e42470aSBarry Smith Input Parameters: 637c7afd0dbSLois Curfman McInnes + snes - the SNES context 638af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6395e42470aSBarry Smith . x - current iterate 6405e42470aSBarry Smith . f - residual evaluated at x 6415e42470aSBarry Smith . y - search direction (contains new iterate on output) 6425e42470aSBarry Smith . w - work vector 643c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6445e42470aSBarry Smith 645c4a48953SLois Curfman McInnes Output Parameters: 646c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6475e42470aSBarry Smith . y - new iterate (contains search direction on input) 6485e42470aSBarry Smith . gnorm - 2-norm of g 6495e42470aSBarry Smith . ynorm - 2-norm of search length 650a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 651fee21e36SBarry Smith 652c4a48953SLois Curfman McInnes Options Database Key: 6534b27c08aSLois Curfman McInnes . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 654c4a48953SLois Curfman McInnes 6555e42470aSBarry Smith Notes: 6564b27c08aSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNESLS method. 6575e42470aSBarry Smith 65836851e7fSLois Curfman McInnes Level: advanced 65936851e7fSLois Curfman McInnes 660f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 661f59ffdedSLois Curfman McInnes 662af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6635e42470aSBarry Smith @*/ 664a7cc72afSBarry 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) 6655e76c431SBarry Smith { 666406556e6SLois Curfman McInnes /* 667406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 668406556e6SLois Curfman McInnes minimization problem: 669406556e6SLois Curfman McInnes min z(x): R^n -> R, 670406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 671406556e6SLois Curfman McInnes */ 6725ca10a37SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 673aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 674ea709b57SSatish Balay PetscScalar cinitslope,clambda; 6756b5873e3SBarry Smith #endif 676dfbe8321SBarry Smith PetscErrorCode ierr; 677a7cc72afSBarry Smith PetscInt count; 6784b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 679ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 6808f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6815e76c431SBarry Smith 6823a40ed3dSBarry Smith PetscFunctionBegin; 683d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 684a7cc72afSBarry Smith *flag = PETSC_TRUE; 6855e42470aSBarry Smith alpha = neP->alpha; 6865e42470aSBarry Smith maxstep = neP->maxstep; 6875e42470aSBarry Smith steptol = neP->steptol; 6885e76c431SBarry Smith 6893505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 69063b9fa5eSBarry Smith if (*ynorm == 0.0) { 691b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 692b37302c6SLois Curfman McInnes *gnorm = fnorm; 693b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 694b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 695a7cc72afSBarry Smith *flag = PETSC_FALSE; 696ad922ac9SBarry Smith goto theend2; 69794a9d846SBarry Smith } 6985e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6995e42470aSBarry Smith scale = maxstep/(*ynorm); 700393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 7015e42470aSBarry Smith *ynorm = maxstep; 7025e76c431SBarry Smith } 703ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 7045ca10a37SBarry Smith minlambda = steptol/rellength; 705a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 706aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 707a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 708329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7095e42470aSBarry Smith #else 710a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7115e42470aSBarry Smith #endif 7125e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7135e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7145e42470aSBarry Smith 715393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 7165334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 71778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 718cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7195d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 720393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 721b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 722ad922ac9SBarry Smith goto theend2; 7235e42470aSBarry Smith } 7245e42470aSBarry Smith 7255e42470aSBarry Smith /* Fit points with quadratic */ 726313b5538SBarry Smith lambda = 1.0; 7275e42470aSBarry Smith count = 1; 7285ca10a37SBarry Smith while (PETSC_TRUE) { 7295e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 73077431f27SBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count); 731b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 732f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 733a7cc72afSBarry Smith *flag = PETSC_FALSE; break; 7345e42470aSBarry Smith } 735a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 736329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 737329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 738329f5518SBarry Smith else lambda = lambdatemp; 739393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7403505fcc1SBarry Smith lambdaneg = -lambda; 741aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7423505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7435e42470aSBarry Smith #else 7443505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7455e42470aSBarry Smith #endif 746*ed98166eSMatthew Knepley if (snes->nfuncs > snes->max_funcs) { 747*ed98166eSMatthew Knepley PetscLogInfo(snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, but unable to find good step length! %D \n",count); 748*ed98166eSMatthew Knepley 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); 749*ed98166eSMatthew Knepley ierr = VecCopy(x,y);CHKERRQ(ierr); 750*ed98166eSMatthew Knepley *flag = PETSC_FALSE; 751*ed98166eSMatthew Knepley break; 752*ed98166eSMatthew Knepley } else { 75378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 754*ed98166eSMatthew Knepley } 755cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7565ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 757393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 758b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 75906259719SBarry Smith break; 7605e42470aSBarry Smith } 7615e42470aSBarry Smith count++; 7625e42470aSBarry Smith } 763ad922ac9SBarry Smith theend2: 7648f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7658f99978dSLois Curfman McInnes if (neP->CheckStep) { 7668f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 767abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 7688f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7698f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7708f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7718f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7728f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7738f99978dSLois Curfman McInnes } 7748f99978dSLois Curfman McInnes } 775d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7763a40ed3dSBarry Smith PetscFunctionReturn(0); 7775e76c431SBarry Smith } 7782343ba6eSBarry Smith 77904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7804a2ae208SSatish Balay #undef __FUNCT__ 7814a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 782c9e6a524SLois Curfman McInnes /*@C 78377c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7844b27c08aSLois Curfman McInnes by the method SNESLS. 7855e76c431SBarry Smith 7865e76c431SBarry Smith Input Parameters: 787c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 788af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 789c7afd0dbSLois Curfman McInnes - func - pointer to int function 7905e76c431SBarry Smith 791fee21e36SBarry Smith Collective on SNES 792fee21e36SBarry Smith 793c4a48953SLois Curfman McInnes Available Routines: 794c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 795f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 796f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 797af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7985e76c431SBarry Smith 799c4a48953SLois Curfman McInnes Options Database Keys: 8004b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8014b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 8024b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 8034b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 8043304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 8053304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 806c4a48953SLois Curfman McInnes 8075e76c431SBarry Smith Calling sequence of func: 808bd208895SLois Curfman McInnes .vb 809af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 810a7cc72afSBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 811bd208895SLois Curfman McInnes .ve 8125e76c431SBarry Smith 8135e76c431SBarry Smith Input parameters for func: 814c7afd0dbSLois Curfman McInnes + snes - nonlinear context 815af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8165e76c431SBarry Smith . x - current iterate 8175e76c431SBarry Smith . f - residual evaluated at x 8185e76c431SBarry Smith . y - search direction (contains new iterate on output) 8195e76c431SBarry Smith . w - work vector 820c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 8215e76c431SBarry Smith 8225e76c431SBarry Smith Output parameters for func: 823c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 8245e76c431SBarry Smith . y - new iterate (contains search direction on input) 8255e76c431SBarry Smith . gnorm - 2-norm of g 8265e76c431SBarry Smith . ynorm - 2-norm of search length 827a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 828f59ffdedSLois Curfman McInnes 82936851e7fSLois Curfman McInnes Level: advanced 83036851e7fSLois Curfman McInnes 831f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 832f59ffdedSLois Curfman McInnes 8335d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 8345d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 8355e76c431SBarry Smith @*/ 836a7cc72afSBarry Smith PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 8375e76c431SBarry Smith { 838a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 83982bf6240SBarry Smith 8403a40ed3dSBarry Smith PetscFunctionBegin; 841c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 84282bf6240SBarry Smith if (f) { 843af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 84482bf6240SBarry Smith } 8453a40ed3dSBarry Smith PetscFunctionReturn(0); 8465e76c431SBarry Smith } 8478e019c35SBarry Smith 848a7cc72afSBarry 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*/ 84904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 850fb2e594dSBarry Smith EXTERN_C_BEGIN 8514a2ae208SSatish Balay #undef __FUNCT__ 8524a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS" 853dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx) 85482bf6240SBarry Smith { 85582bf6240SBarry Smith PetscFunctionBegin; 8564b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 8574b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 85882bf6240SBarry Smith PetscFunctionReturn(0); 85982bf6240SBarry Smith } 860fb2e594dSBarry Smith EXTERN_C_END 86104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8624a2ae208SSatish Balay #undef __FUNCT__ 8634a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck" 864c8dd0c0dSLois Curfman McInnes /*@C 865530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 8664b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 867c8dd0c0dSLois Curfman McInnes 868c8dd0c0dSLois Curfman McInnes Input Parameters: 869c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 870c8dd0c0dSLois Curfman McInnes . func - pointer to int function 871c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 872c8dd0c0dSLois Curfman McInnes 873c8dd0c0dSLois Curfman McInnes Collective on SNES 874c8dd0c0dSLois Curfman McInnes 875c8dd0c0dSLois Curfman McInnes Calling sequence of func: 876c8dd0c0dSLois Curfman McInnes .vb 877b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 878c8dd0c0dSLois Curfman McInnes .ve 879b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 880b1ae27eaSLois Curfman McInnes on failure. 881c8dd0c0dSLois Curfman McInnes 882c8dd0c0dSLois Curfman McInnes Input parameters for func: 883c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 884c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8858f99978dSLois Curfman McInnes - x - current candidate iterate 886c8dd0c0dSLois Curfman McInnes 887c8dd0c0dSLois Curfman McInnes Output parameters for func: 888c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 889c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 890c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 891c8dd0c0dSLois Curfman McInnes 892c8dd0c0dSLois Curfman McInnes Level: advanced 893c8dd0c0dSLois Curfman McInnes 8948f99978dSLois Curfman McInnes Notes: 895b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 896b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 897b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 898b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 899ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9008f99978dSLois Curfman McInnes 901b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 902b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 903b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 904ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 905ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 906ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 907ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 908b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9098f99978dSLois Curfman McInnes 910c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 911c8dd0c0dSLois Curfman McInnes 912c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 913c8dd0c0dSLois Curfman McInnes @*/ 9146849ba73SBarry Smith PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 915c8dd0c0dSLois Curfman McInnes { 9166849ba73SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*); 917c8dd0c0dSLois Curfman McInnes 918c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 919c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 920c8dd0c0dSLois Curfman McInnes if (f) { 921c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 922c8dd0c0dSLois Curfman McInnes } 923c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 924c8dd0c0dSLois Curfman McInnes } 925c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9266849ba73SBarry Smith typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/ 927c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 9284a2ae208SSatish Balay #undef __FUNCT__ 9294a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS" 930dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx) 931c8dd0c0dSLois Curfman McInnes { 932c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 9334b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 9344b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 935c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 936c8dd0c0dSLois Curfman McInnes } 937c8dd0c0dSLois Curfman McInnes EXTERN_C_END 938c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 93904d965bbSLois Curfman McInnes /* 9404b27c08aSLois Curfman McInnes SNESPrintHelp_LS - Prints all options for the SNES_LS method. 941329e7e40SMatthew Knepley 942329e7e40SMatthew Knepley Input Parameter: 943329e7e40SMatthew Knepley . snes - the SNES context 944329e7e40SMatthew Knepley 945329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 946329e7e40SMatthew Knepley */ 947329e7e40SMatthew Knepley #undef __FUNCT__ 9484b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS" 9496849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p) 950329e7e40SMatthew Knepley { 9514b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 952329e7e40SMatthew Knepley 953329e7e40SMatthew Knepley PetscFunctionBegin; 9544b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 9554b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 9564b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 9574b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 9584b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 959329e7e40SMatthew Knepley PetscFunctionReturn(0); 960329e7e40SMatthew Knepley } 961329e7e40SMatthew Knepley 962329e7e40SMatthew Knepley /* 9634b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 96404d965bbSLois Curfman McInnes 96504d965bbSLois Curfman McInnes Input Parameters: 96604d965bbSLois Curfman McInnes . SNES - the SNES context 96704d965bbSLois Curfman McInnes . viewer - visualization context 96804d965bbSLois Curfman McInnes 96904d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 97004d965bbSLois Curfman McInnes */ 9714a2ae208SSatish Balay #undef __FUNCT__ 9724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 9736849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 974a935fc98SLois Curfman McInnes { 9754b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 9762fc52814SBarry Smith const char *cstr; 977dfbe8321SBarry Smith PetscErrorCode ierr; 97832077d6dSBarry Smith PetscTruth iascii; 979a935fc98SLois Curfman McInnes 9803a40ed3dSBarry Smith PetscFunctionBegin; 98132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 98232077d6dSBarry Smith if (iascii) { 98319bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 98419bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 98519bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 98619bcc07fSBarry Smith else cstr = "unknown"; 987b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 988b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9895cd90555SBarry Smith } else { 99079a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 99119bcc07fSBarry Smith } 9923a40ed3dSBarry Smith PetscFunctionReturn(0); 993a935fc98SLois Curfman McInnes } 99404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 99504d965bbSLois Curfman McInnes /* 9964b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 99704d965bbSLois Curfman McInnes 99804d965bbSLois Curfman McInnes Input Parameter: 99904d965bbSLois Curfman McInnes . snes - the SNES context 100004d965bbSLois Curfman McInnes 100104d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 100204d965bbSLois Curfman McInnes */ 10034a2ae208SSatish Balay #undef __FUNCT__ 10044b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 10056849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 10065e42470aSBarry Smith { 10074b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1008e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1009dfbe8321SBarry Smith PetscErrorCode ierr; 1010a7cc72afSBarry Smith PetscInt indx; 1011f1af5d2fSBarry Smith PetscTruth flg; 10125e42470aSBarry Smith 10133a40ed3dSBarry Smith PetscFunctionBegin; 1014b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 10154b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 10164b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 10174b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1018186905e3SBarry Smith 101922e36657SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 102025cdf11fSBarry Smith if (flg) { 102122e36657SBarry Smith switch (indx) { 1022b49fd9e1SBarry Smith case 0: 1023af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1024b49fd9e1SBarry Smith break; 1025b49fd9e1SBarry Smith case 1: 1026af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1027b49fd9e1SBarry Smith break; 1028b49fd9e1SBarry Smith case 2: 1029af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1030b49fd9e1SBarry Smith break; 1031b49fd9e1SBarry Smith case 3: 1032af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1033b49fd9e1SBarry Smith break; 10345e42470aSBarry Smith } 10355e42470aSBarry Smith } 1036b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10373a40ed3dSBarry Smith PetscFunctionReturn(0); 10385e42470aSBarry Smith } 103904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1040ebe3b25bSBarry Smith /*MC 1041ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 104204d965bbSLois Curfman McInnes 1043ebe3b25bSBarry Smith Options Database: 1044ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1045ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1046ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1047ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1048ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1049ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 105004d965bbSLois Curfman McInnes 1051ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 105204d965bbSLois Curfman McInnes 1053ee3001cbSBarry Smith Level: beginner 1054ee3001cbSBarry Smith 1055ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1056ebe3b25bSBarry Smith SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1057ebe3b25bSBarry Smith SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1058ebe3b25bSBarry Smith 1059ebe3b25bSBarry Smith M*/ 1060fb2e594dSBarry Smith EXTERN_C_BEGIN 10614a2ae208SSatish Balay #undef __FUNCT__ 10624b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 1063dfbe8321SBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 10645e42470aSBarry Smith { 1065dfbe8321SBarry Smith PetscErrorCode ierr; 10664b27c08aSLois Curfman McInnes SNES_LS *neP; 10675e42470aSBarry Smith 10683a40ed3dSBarry Smith PetscFunctionBegin; 10694b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 10704b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 10714b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 10724b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 10734b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 10744b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 10754b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 10765baf8537SBarry Smith snes->nwork = 0; 10775e42470aSBarry Smith 10784b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 107952e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 10805e42470aSBarry Smith snes->data = (void*)neP; 10815e42470aSBarry Smith neP->alpha = 1.e-4; 10825e42470aSBarry Smith neP->maxstep = 1.e8; 10835e42470aSBarry Smith neP->steptol = 1.e-12; 10845e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1085c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1086c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1087c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 108882bf6240SBarry Smith 10895d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10905d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 109182bf6240SBarry Smith 10923a40ed3dSBarry Smith PetscFunctionReturn(0); 10935e42470aSBarry Smith } 1094fb2e594dSBarry Smith EXTERN_C_END 109582bf6240SBarry Smith 109682bf6240SBarry Smith 109782bf6240SBarry Smith 1098