173f4d377SMatthew Knepley /*$Id: ls.c,v 1.172 2001/08/07 03:04:11 balay Exp $*/ 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 58a5d9ee4SBarry Smith /* 68a5d9ee4SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the function, 78a5d9ee4SBarry Smith but not a zero. In the case when one cannot compute J^T F we use the fact that 836669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 936669109SBarry Smith for this trick. 108a5d9ee4SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 138a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 168a5d9ee4SBarry Smith int ierr; 1736669109SBarry Smith PetscTruth hastranspose; 188a5d9ee4SBarry Smith 198a5d9ee4SBarry Smith PetscFunctionBegin; 208a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2136669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2236669109SBarry Smith if (hastranspose) { 238a5d9ee4SBarry Smith /* Compute || J^T F|| */ 248a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 258a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 264b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 30ea709b57SSatish Balay PetscScalar result; 3136669109SBarry Smith PetscReal wnorm; 3236669109SBarry Smith 3336669109SBarry Smith ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDestroy(work);CHKERRQ(ierr); 3936669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 404b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 475ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4874637425SBarry Smith */ 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 5174637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5274637425SBarry Smith { 5374637425SBarry Smith PetscReal a1,a2; 5474637425SBarry Smith int ierr; 5574637425SBarry Smith PetscTruth hastranspose; 56ea709b57SSatish Balay PetscScalar mone = -1.0; 5774637425SBarry Smith 5874637425SBarry Smith PetscFunctionBegin; 5974637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 6074637425SBarry Smith if (hastranspose) { 6174637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6274637425SBarry Smith ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr); 6374637425SBarry Smith 6474637425SBarry Smith /* Compute || J^T W|| */ 6574637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6774637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 689994e62eSSatish Balay if (a1 != 0) { 694b27c08aSLois Curfman McInnes PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1); 7085471664SBarry Smith } 7174637425SBarry Smith } 7274637425SBarry Smith PetscFunctionReturn(0); 7374637425SBarry Smith } 7474637425SBarry Smith 7504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7604d965bbSLois Curfman McInnes 7704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7894b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8004d965bbSLois Curfman McInnes respectively. 8104d965bbSLois Curfman McInnes 82fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 85fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 884b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 894b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9304d965bbSLois Curfman McInnes for all nonlinear solvers. 9404d965bbSLois Curfman McInnes 9504d965bbSLois Curfman McInnes Another key routine is: 9604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9904d965bbSLois Curfman McInnes SNESSolve() if necessary. 10004d965bbSLois Curfman McInnes 10104d965bbSLois Curfman McInnes Additional basic routines are: 10204d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10304d965bbSLois Curfman McInnes have actually been used. 10404d965bbSLois Curfman McInnes These are called by application codes via the interface routines 105186905e3SBarry Smith SNESView(). 10604d965bbSLois Curfman McInnes 10704d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10804d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10904d965bbSLois Curfman McInnes above description applies to these categories also. 11004d965bbSLois Curfman McInnes 11104d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1125e42470aSBarry Smith /* 1134b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11404d965bbSLois Curfman McInnes method with a line search. 1155e76c431SBarry Smith 11604d965bbSLois Curfman McInnes Input Parameters: 11704d965bbSLois Curfman McInnes . snes - the SNES context 1185e76c431SBarry Smith 11904d965bbSLois Curfman McInnes Output Parameter: 120c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12104d965bbSLois Curfman McInnes 12204d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1235e76c431SBarry Smith 1245e76c431SBarry Smith Notes: 1255e76c431SBarry Smith This implements essentially a truncated Newton method with a 1265e76c431SBarry Smith line search. By default a cubic backtracking line search 1275e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1285e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 129393d2d9aSLois Curfman McInnes and Schnabel. 1305e76c431SBarry Smith */ 1314a2ae208SSatish Balay #undef __FUNCT__ 1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 133c9780b6fSBarry Smith int SNESSolve_LS(SNES snes) 1345e76c431SBarry Smith { 1354b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 13684c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 137112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 138329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1395e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 140c293cc10SBarry Smith KSP ksp; 1415e76c431SBarry Smith 1423a40ed3dSBarry Smith PetscFunctionBegin; 14394b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 144184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 145184914b5SBarry Smith 1465e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1475e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1495e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1505e42470aSBarry Smith G = snes->work[1]; 1515e42470aSBarry Smith W = snes->work[2]; 1525e76c431SBarry Smith 1534c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1544c49b128SBarry Smith snes->iter = 0; 1554c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1565334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 157cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1580f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1595e42470aSBarry Smith snes->norm = fnorm; 1600f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16184c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16294a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1635e76c431SBarry Smith 164c9780b6fSBarry Smith if (fnorm < snes->atol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16594a9d846SBarry Smith 166d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 167d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 168d034289bSBarry Smith 1695e76c431SBarry Smith for (i=0; i<maxits; i++) { 1705e76c431SBarry Smith 171329e7e40SMatthew Knepley /* Call general purpose update function */ 172de8cb200SMatthew Knepley if (snes->update != PETSC_NULL) { 173329e7e40SMatthew Knepley ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr); 174de8cb200SMatthew Knepley } 175329e7e40SMatthew Knepley 176ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1775334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 17894b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 17994b7f48cSBarry Smith ierr = KSPSetRhs(snes->ksp,F);CHKERRQ(ierr); 18094b7f48cSBarry Smith ierr = KSPSetSolution(snes->ksp,Y);CHKERRQ(ierr); 18194b7f48cSBarry Smith ierr = KSPSolve(snes->ksp);CHKERRQ(ierr); 182c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 18374637425SBarry Smith 184b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 18574637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 18685471664SBarry Smith } 18774637425SBarry Smith 18890cbd584SBarry Smith /* should check what happened to the linear solve? */ 1893505fcc1SBarry Smith snes->linear_its += lits; 1904b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 191ea4d3ed3SLois Curfman McInnes 192ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 193ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 194ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 195ea4d3ed3SLois Curfman McInnes */ 19681b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 197af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 1984b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 199a3b891d8SBarry Smith 200a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 201a3b891d8SBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 202a3b891d8SBarry Smith fnorm = gnorm; 203a3b891d8SBarry Smith 2045ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2055ed2d596SBarry Smith snes->iter = i+1; 2065ed2d596SBarry Smith snes->norm = fnorm; 2075ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2085ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2095ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2105ed2d596SBarry Smith 2113505fcc1SBarry Smith if (lsfail) { 2128a5d9ee4SBarry Smith PetscTruth ismin; 21350ffb88aSMatthew Knepley 21450ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2153505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2168a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2178a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2183505fcc1SBarry Smith break; 2193505fcc1SBarry Smith } 22050ffb88aSMatthew Knepley } 2215e76c431SBarry Smith 2225e76c431SBarry Smith /* Test for convergence */ 22329e0b56fSBarry Smith if (snes->converged) { 22429e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2253505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2263505fcc1SBarry Smith if (snes->reason) { 22716c95cacSBarry Smith break; 22816c95cacSBarry Smith } 22916c95cacSBarry Smith } 23029e0b56fSBarry Smith } 23139e2f89bSBarry Smith if (X != snes->vec_sol) { 232393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 233e15f7bb5SBarry Smith } 234e15f7bb5SBarry Smith if (F != snes->vec_func) { 235e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 236e15f7bb5SBarry Smith } 23739e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 23839e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 23952392280SLois Curfman McInnes if (i == maxits) { 2404b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %d\n",maxits); 2413505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 24252392280SLois Curfman McInnes } 2433a40ed3dSBarry Smith PetscFunctionReturn(0); 2445e76c431SBarry Smith } 24504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 24604d965bbSLois Curfman McInnes /* 2474b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2484b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 24904d965bbSLois Curfman McInnes 25004d965bbSLois Curfman McInnes Input Parameter: 25104d965bbSLois Curfman McInnes . snes - the SNES context 25204d965bbSLois Curfman McInnes . x - the solution vector 25304d965bbSLois Curfman McInnes 25404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 25504d965bbSLois Curfman McInnes 25604d965bbSLois Curfman McInnes Notes: 2574b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 25804d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 25904d965bbSLois Curfman McInnes the call to SNESSolve(). 26004d965bbSLois Curfman McInnes */ 2614a2ae208SSatish Balay #undef __FUNCT__ 2624b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 2634b27c08aSLois Curfman McInnes int SNESSetUp_LS(SNES snes) 2645e76c431SBarry Smith { 2655e42470aSBarry Smith int ierr; 2663a40ed3dSBarry Smith 2673a40ed3dSBarry Smith PetscFunctionBegin; 26881b6cf68SLois Curfman McInnes snes->nwork = 4; 269d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 270b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 27181b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2723a40ed3dSBarry Smith PetscFunctionReturn(0); 2735e76c431SBarry Smith } 27404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27504d965bbSLois Curfman McInnes /* 2764b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2774b27c08aSLois Curfman McInnes with SNESCreate_LS(). 27804d965bbSLois Curfman McInnes 27904d965bbSLois Curfman McInnes Input Parameter: 28004d965bbSLois Curfman McInnes . snes - the SNES context 28104d965bbSLois Curfman McInnes 282de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 28304d965bbSLois Curfman McInnes */ 2844a2ae208SSatish Balay #undef __FUNCT__ 2854b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 2864b27c08aSLois Curfman McInnes int SNESDestroy_LS(SNES snes) 2875e76c431SBarry Smith { 288393d2d9aSLois Curfman McInnes int ierr; 2893a40ed3dSBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 2915baf8537SBarry Smith if (snes->nwork) { 2924b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 29321c89e3eSBarry Smith } 2945c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2953a40ed3dSBarry Smith PetscFunctionReturn(0); 2965e76c431SBarry Smith } 29704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2984a2ae208SSatish Balay #undef __FUNCT__ 2994a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch" 30082bf6240SBarry Smith 3014b828684SBarry Smith /*@C 3025e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 3035e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3045e76c431SBarry Smith to serve as a template and is not recommended for general use. 3055e76c431SBarry Smith 306c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 307c7afd0dbSLois Curfman McInnes 3085e76c431SBarry Smith Input Parameters: 309c7afd0dbSLois Curfman McInnes + snes - nonlinear context 310af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3115e76c431SBarry Smith . x - current iterate 3125e76c431SBarry Smith . f - residual evaluated at x 3135e76c431SBarry Smith . y - search direction (contains new iterate on output) 3145e76c431SBarry Smith . w - work vector 315c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3165e76c431SBarry Smith 317c4a48953SLois Curfman McInnes Output Parameters: 318c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3195e76c431SBarry Smith . y - new iterate (contains search direction on input) 3205e76c431SBarry Smith . gnorm - 2-norm of g 3215e76c431SBarry Smith . ynorm - 2-norm of search length 322c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 323fee21e36SBarry Smith 324c4a48953SLois Curfman McInnes Options Database Key: 3254b27c08aSLois Curfman McInnes . -snes_ls basic - Activates SNESNoLineSearch() 326c4a48953SLois Curfman McInnes 32736851e7fSLois Curfman McInnes Level: advanced 32836851e7fSLois Curfman McInnes 32928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 33028ae5a21SLois Curfman McInnes 331f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 332af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3335e76c431SBarry Smith @*/ 3345d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 3355e76c431SBarry Smith { 3365e42470aSBarry Smith int ierr; 337ea709b57SSatish Balay PetscScalar mone = -1.0; 3384b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3398f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3405334005bSBarry Smith 3413a40ed3dSBarry Smith PetscFunctionBegin; 342761aaf1bSLois Curfman McInnes *flag = 0; 343d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 344a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 345ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3468f99978dSLois Curfman McInnes if (neP->CheckStep) { 3478f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3488f99978dSLois Curfman McInnes } 349ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 350a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 351d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3523a40ed3dSBarry Smith PetscFunctionReturn(0); 3535e76c431SBarry Smith } 35404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3554a2ae208SSatish Balay #undef __FUNCT__ 3564a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms" 35782bf6240SBarry Smith 35829e0b56fSBarry Smith /*@C 35929e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 36029e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 36129e0b56fSBarry Smith even compute the norm of the function or search direction; this 36229e0b56fSBarry Smith is intended only when you know the full step is fine and are 36329e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 364c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 365c7afd0dbSLois Curfman McInnes 366c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 36729e0b56fSBarry Smith 36829e0b56fSBarry Smith Input Parameters: 369c7afd0dbSLois Curfman McInnes + snes - nonlinear context 370af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 37129e0b56fSBarry Smith . x - current iterate 37229e0b56fSBarry Smith . f - residual evaluated at x 37329e0b56fSBarry Smith . y - search direction (contains new iterate on output) 37429e0b56fSBarry Smith . w - work vector 375c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 37629e0b56fSBarry Smith 37729e0b56fSBarry Smith Output Parameters: 378c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 37929e0b56fSBarry Smith . gnorm - not changed 38029e0b56fSBarry Smith . ynorm - not changed 381c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 382fee21e36SBarry Smith 38329e0b56fSBarry Smith Options Database Key: 3844b27c08aSLois Curfman McInnes . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 38529e0b56fSBarry Smith 3868cbba510SBarry Smith Notes: 387ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 388ea56f5baSLois Curfman McInnes either the options 389ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 390ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 391329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 392329f5518SBarry Smith otherwise, the SNES solver will generate an error. 393329f5518SBarry Smith 394329f5518SBarry Smith During the final iteration this will not evaluate the function at 395329f5518SBarry Smith the solution point. This is to save a function evaluation while 396329f5518SBarry Smith using pseudo-timestepping. 3978cbba510SBarry Smith 398ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 399ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 400ea56f5baSLois Curfman McInnes correct, since they are not computed. 401ea56f5baSLois Curfman McInnes 402ea56f5baSLois Curfman McInnes Level: advanced 4038cbba510SBarry Smith 40429e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 40529e0b56fSBarry Smith 40629e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 40729e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 40829e0b56fSBarry Smith @*/ 4095d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 41029e0b56fSBarry Smith { 41129e0b56fSBarry Smith int ierr; 412ea709b57SSatish Balay PetscScalar mone = -1.0; 4134b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4148f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 41529e0b56fSBarry Smith 4163a40ed3dSBarry Smith PetscFunctionBegin; 4178cbba510SBarry Smith *flag = 0; 418d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 41929e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4208f99978dSLois Curfman McInnes if (neP->CheckStep) { 4218f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4228f99978dSLois Curfman McInnes } 423329f5518SBarry Smith 424329f5518SBarry Smith /* don't evaluate function the last time through */ 425329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 42629e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 427329f5518SBarry Smith } 428d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4293a40ed3dSBarry Smith PetscFunctionReturn(0); 43029e0b56fSBarry Smith } 43104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4324a2ae208SSatish Balay #undef __FUNCT__ 4334a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch" 4344b828684SBarry Smith /*@C 435f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4365e76c431SBarry Smith 437c7afd0dbSLois Curfman McInnes Collective on SNES 438c7afd0dbSLois Curfman McInnes 4395e76c431SBarry Smith Input Parameters: 440c7afd0dbSLois Curfman McInnes + snes - nonlinear context 441af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4425e76c431SBarry Smith . x - current iterate 4435e76c431SBarry Smith . f - residual evaluated at x 4445e76c431SBarry Smith . y - search direction (contains new iterate on output) 4455e76c431SBarry Smith . w - work vector 446c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4475e76c431SBarry Smith 448393d2d9aSLois Curfman McInnes Output Parameters: 449c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4505e76c431SBarry Smith . y - new iterate (contains search direction on input) 4515e76c431SBarry Smith . gnorm - 2-norm of g 4525e76c431SBarry Smith . ynorm - 2-norm of search length 453c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 454fee21e36SBarry Smith 455c4a48953SLois Curfman McInnes Options Database Key: 4564b27c08aSLois Curfman McInnes $ -snes_ls cubic - Activates SNESCubicLineSearch() 457c4a48953SLois Curfman McInnes 4585e76c431SBarry Smith Notes: 4595e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4605e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4615e76c431SBarry Smith 46236851e7fSLois Curfman McInnes Level: advanced 46336851e7fSLois Curfman McInnes 46428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 46528ae5a21SLois Curfman McInnes 466af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4675e76c431SBarry Smith @*/ 4685d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 4695e76c431SBarry Smith { 470406556e6SLois Curfman McInnes /* 471406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 472406556e6SLois Curfman McInnes minimization problem: 473406556e6SLois Curfman McInnes min z(x): R^n -> R, 474406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 475406556e6SLois Curfman McInnes */ 476406556e6SLois Curfman McInnes 4775ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 478329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 479aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 480ea709b57SSatish Balay PetscScalar cinitslope,clambda; 4816b5873e3SBarry Smith #endif 4825e42470aSBarry Smith int ierr,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); 489761aaf1bSLois Curfman McInnes *flag = 0; 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); 50063b9fa5eSBarry Smith *flag = -1; 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 */ 562b0a32e0cSBarry 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); 565761aaf1bSLois Curfman McInnes *flag = -1; 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); 6098f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* 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 642c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 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 @*/ 6565d1a10b1SBarry Smith int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 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 6685e42470aSBarry Smith int ierr,count; 6694b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 670ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 6718f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6725e76c431SBarry Smith 6733a40ed3dSBarry Smith PetscFunctionBegin; 674d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 675761aaf1bSLois Curfman McInnes *flag = 0; 6765e42470aSBarry Smith alpha = neP->alpha; 6775e42470aSBarry Smith maxstep = neP->maxstep; 6785e42470aSBarry Smith steptol = neP->steptol; 6795e76c431SBarry Smith 6803505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 68163b9fa5eSBarry Smith if (*ynorm == 0.0) { 682b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 683b37302c6SLois Curfman McInnes *gnorm = fnorm; 684b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 685b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 68663b9fa5eSBarry Smith *flag = -1; 687ad922ac9SBarry Smith goto theend2; 68894a9d846SBarry Smith } 6895e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6905e42470aSBarry Smith scale = maxstep/(*ynorm); 691393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6925e42470aSBarry Smith *ynorm = maxstep; 6935e76c431SBarry Smith } 694ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 6955ca10a37SBarry Smith minlambda = steptol/rellength; 696a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 697aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 698a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 699329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7005e42470aSBarry Smith #else 701a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7025e42470aSBarry Smith #endif 7035e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7045e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7055e42470aSBarry Smith 706393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 7075334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 70878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 709cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7105d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 711393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 712b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 713ad922ac9SBarry Smith goto theend2; 7145e42470aSBarry Smith } 7155e42470aSBarry Smith 7165e42470aSBarry Smith /* Fit points with quadratic */ 717313b5538SBarry Smith lambda = 1.0; 7185e42470aSBarry Smith count = 1; 7195ca10a37SBarry Smith while (PETSC_TRUE) { 7205e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 721b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 722b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 723f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 724761aaf1bSLois Curfman McInnes *flag = -1; break; 7255e42470aSBarry Smith } 726a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 727329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 728329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 729329f5518SBarry Smith else lambda = lambdatemp; 730393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7313505fcc1SBarry Smith lambdaneg = -lambda; 732aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7333505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7345e42470aSBarry Smith #else 7353505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7365e42470aSBarry Smith #endif 73778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 738cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7395ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 740393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 741b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 74206259719SBarry Smith break; 7435e42470aSBarry Smith } 7445e42470aSBarry Smith count++; 7455e42470aSBarry Smith } 746ad922ac9SBarry Smith theend2: 7478f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7488f99978dSLois Curfman McInnes if (neP->CheckStep) { 7498f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7508f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7518f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7528f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7538f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7548f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7558f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7568f99978dSLois Curfman McInnes } 7578f99978dSLois Curfman McInnes } 758d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7593a40ed3dSBarry Smith PetscFunctionReturn(0); 7605e76c431SBarry Smith } 7612343ba6eSBarry Smith 76204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7634a2ae208SSatish Balay #undef __FUNCT__ 7644a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 765c9e6a524SLois Curfman McInnes /*@C 76677c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7674b27c08aSLois Curfman McInnes by the method SNESLS. 7685e76c431SBarry Smith 7695e76c431SBarry Smith Input Parameters: 770c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 771af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 772c7afd0dbSLois Curfman McInnes - func - pointer to int function 7735e76c431SBarry Smith 774fee21e36SBarry Smith Collective on SNES 775fee21e36SBarry Smith 776c4a48953SLois Curfman McInnes Available Routines: 777c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 778f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 779f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 780af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7815e76c431SBarry Smith 782c4a48953SLois Curfman McInnes Options Database Keys: 7834b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 7844b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 7854b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 7864b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 7873304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 7883304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 789c4a48953SLois Curfman McInnes 7905e76c431SBarry Smith Calling sequence of func: 791bd208895SLois Curfman McInnes .vb 792af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 793329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 794bd208895SLois Curfman McInnes .ve 7955e76c431SBarry Smith 7965e76c431SBarry Smith Input parameters for func: 797c7afd0dbSLois Curfman McInnes + snes - nonlinear context 798af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7995e76c431SBarry Smith . x - current iterate 8005e76c431SBarry Smith . f - residual evaluated at x 8015e76c431SBarry Smith . y - search direction (contains new iterate on output) 8025e76c431SBarry Smith . w - work vector 803c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 8045e76c431SBarry Smith 8055e76c431SBarry Smith Output parameters for func: 806c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 8075e76c431SBarry Smith . y - new iterate (contains search direction on input) 8085e76c431SBarry Smith . gnorm - 2-norm of g 8095e76c431SBarry Smith . ynorm - 2-norm of search length 810c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 811761aaf1bSLois Curfman McInnes 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 @*/ 820329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 8215e76c431SBarry Smith { 822329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),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 8328e019c35SBarry Smith typedef int (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*); /* 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" 8378e019c35SBarry Smith int 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 @*/ 8988f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 899c8dd0c0dSLois Curfman McInnes { 9008f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(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 /* -------------------------------------------------------------------------- */ 9102343ba6eSBarry Smith typedef int (*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" 9142343ba6eSBarry Smith int 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" 9334b27c08aSLois Curfman McInnes static int 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" 9574b27c08aSLois Curfman McInnes static int SNESView_LS(SNES snes,PetscViewer viewer) 958a935fc98SLois Curfman McInnes { 9594b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 9602fc52814SBarry Smith const char *cstr; 96151695f54SLois Curfman McInnes int ierr; 9626831982aSBarry Smith PetscTruth isascii; 963a935fc98SLois Curfman McInnes 9643a40ed3dSBarry Smith PetscFunctionBegin; 965b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 9660f5bd95cSBarry Smith if (isascii) { 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 { 97429bbc08cSBarry Smith SETERRQ1(1,"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" 9894b27c08aSLois Curfman McInnes static int 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"}; 99322e36657SBarry Smith int ierr,indx; 994f1af5d2fSBarry Smith PetscTruth flg; 9955e42470aSBarry Smith 9963a40ed3dSBarry Smith PetscFunctionBegin; 997b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 9984b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 9994b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 10004b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1001186905e3SBarry Smith 100222e36657SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 100325cdf11fSBarry Smith if (flg) { 100422e36657SBarry Smith switch (indx) { 1005b49fd9e1SBarry Smith case 0: 1006af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1007b49fd9e1SBarry Smith break; 1008b49fd9e1SBarry Smith case 1: 1009af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1010b49fd9e1SBarry Smith break; 1011b49fd9e1SBarry Smith case 2: 1012af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1013b49fd9e1SBarry Smith break; 1014b49fd9e1SBarry Smith case 3: 1015af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1016b49fd9e1SBarry Smith break; 10175e42470aSBarry Smith } 10185e42470aSBarry Smith } 1019b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10203a40ed3dSBarry Smith PetscFunctionReturn(0); 10215e42470aSBarry Smith } 102204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1023ebe3b25bSBarry Smith /*MC 1024ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 102504d965bbSLois Curfman McInnes 1026ebe3b25bSBarry Smith Options Database: 1027ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1028ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1029ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1030ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1031ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1032ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 103304d965bbSLois Curfman McInnes 1034ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 103504d965bbSLois Curfman McInnes 1036*ee3001cbSBarry Smith Level: beginner 1037*ee3001cbSBarry Smith 1038ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1039ebe3b25bSBarry Smith SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1040ebe3b25bSBarry Smith SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1041ebe3b25bSBarry Smith 1042ebe3b25bSBarry Smith M*/ 1043fb2e594dSBarry Smith EXTERN_C_BEGIN 10444a2ae208SSatish Balay #undef __FUNCT__ 10454b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 10464b27c08aSLois Curfman McInnes int SNESCreate_LS(SNES snes) 10475e42470aSBarry Smith { 104882bf6240SBarry Smith int ierr; 10494b27c08aSLois Curfman McInnes SNES_LS *neP; 10505e42470aSBarry Smith 10513a40ed3dSBarry Smith PetscFunctionBegin; 10524b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 10534b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 10544b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 10554b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 10564b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 10574b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 10584b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 10595baf8537SBarry Smith snes->nwork = 0; 10605e42470aSBarry Smith 10614b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 10624b27c08aSLois Curfman McInnes PetscLogObjectMemory(snes,sizeof(SNES_LS)); 10635e42470aSBarry Smith snes->data = (void*)neP; 10645e42470aSBarry Smith neP->alpha = 1.e-4; 10655e42470aSBarry Smith neP->maxstep = 1.e8; 10665e42470aSBarry Smith neP->steptol = 1.e-12; 10675e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1068c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1069c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1070c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 107182bf6240SBarry Smith 10725d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10735d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 107482bf6240SBarry Smith 10753a40ed3dSBarry Smith PetscFunctionReturn(0); 10765e42470aSBarry Smith } 1077fb2e594dSBarry Smith EXTERN_C_END 107882bf6240SBarry Smith 107982bf6240SBarry Smith 108082bf6240SBarry Smith 1081