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, 7804d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 7904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8004d965bbSLois Curfman McInnes respectively. 8104d965bbSLois Curfman McInnes 82fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 85fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 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" 1334b27c08aSLois Curfman McInnes int SNESSolve_LS(SNES snes,int *outits) 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; 140*c293cc10SBarry Smith KSP ksp; 1415e76c431SBarry Smith 1423a40ed3dSBarry Smith PetscFunctionBegin; 143*c293cc10SBarry Smith ierr = SLESGetKSP(snes->sles,&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 164184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; 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); 1785334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 179*c293cc10SBarry Smith ierr = SLESSolve(snes->sles,F,Y);CHKERRQ(ierr); 180*c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 18174637425SBarry Smith 182b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 18374637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 18485471664SBarry Smith } 18574637425SBarry Smith 18690cbd584SBarry Smith /* should check what happened to the linear solve? */ 1873505fcc1SBarry Smith snes->linear_its += lits; 1884b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 189ea4d3ed3SLois Curfman McInnes 190ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 191ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 192ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 193ea4d3ed3SLois Curfman McInnes */ 19481b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 195af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 1964b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 197a3b891d8SBarry Smith 198a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 199a3b891d8SBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 200a3b891d8SBarry Smith fnorm = gnorm; 201a3b891d8SBarry Smith 2025ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2035ed2d596SBarry Smith snes->iter = i+1; 2045ed2d596SBarry Smith snes->norm = fnorm; 2055ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2065ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2075ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2085ed2d596SBarry Smith 2093505fcc1SBarry Smith if (lsfail) { 2108a5d9ee4SBarry Smith PetscTruth ismin; 21150ffb88aSMatthew Knepley 21250ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2133505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2148a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2158a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2163505fcc1SBarry Smith break; 2173505fcc1SBarry Smith } 21850ffb88aSMatthew Knepley } 2195e76c431SBarry Smith 2205e76c431SBarry Smith /* Test for convergence */ 22129e0b56fSBarry Smith if (snes->converged) { 22229e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2233505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2243505fcc1SBarry Smith if (snes->reason) { 22516c95cacSBarry Smith break; 22616c95cacSBarry Smith } 22716c95cacSBarry Smith } 22829e0b56fSBarry Smith } 22939e2f89bSBarry Smith if (X != snes->vec_sol) { 230393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 231e15f7bb5SBarry Smith } 232e15f7bb5SBarry Smith if (F != snes->vec_func) { 233e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 234e15f7bb5SBarry Smith } 23539e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 23639e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 23752392280SLois Curfman McInnes if (i == maxits) { 2384b27c08aSLois Curfman McInnes PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %d\n",maxits); 23952392280SLois Curfman McInnes i--; 2403505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 24152392280SLois Curfman McInnes } 2420f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2430f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2445e42470aSBarry Smith *outits = i+1; 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2465e76c431SBarry Smith } 24704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 24804d965bbSLois Curfman McInnes /* 2494b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2504b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 25104d965bbSLois Curfman McInnes 25204d965bbSLois Curfman McInnes Input Parameter: 25304d965bbSLois Curfman McInnes . snes - the SNES context 25404d965bbSLois Curfman McInnes . x - the solution vector 25504d965bbSLois Curfman McInnes 25604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 25704d965bbSLois Curfman McInnes 25804d965bbSLois Curfman McInnes Notes: 2594b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 26004d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 26104d965bbSLois Curfman McInnes the call to SNESSolve(). 26204d965bbSLois Curfman McInnes */ 2634a2ae208SSatish Balay #undef __FUNCT__ 2644b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 2654b27c08aSLois Curfman McInnes int SNESSetUp_LS(SNES snes) 2665e76c431SBarry Smith { 2675e42470aSBarry Smith int ierr; 2683a40ed3dSBarry Smith 2693a40ed3dSBarry Smith PetscFunctionBegin; 27081b6cf68SLois Curfman McInnes snes->nwork = 4; 271d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 272b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 27381b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2743a40ed3dSBarry Smith PetscFunctionReturn(0); 2755e76c431SBarry Smith } 27604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27704d965bbSLois Curfman McInnes /* 2784b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2794b27c08aSLois Curfman McInnes with SNESCreate_LS(). 28004d965bbSLois Curfman McInnes 28104d965bbSLois Curfman McInnes Input Parameter: 28204d965bbSLois Curfman McInnes . snes - the SNES context 28304d965bbSLois Curfman McInnes 284de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 28504d965bbSLois Curfman McInnes */ 2864a2ae208SSatish Balay #undef __FUNCT__ 2874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 2884b27c08aSLois Curfman McInnes int SNESDestroy_LS(SNES snes) 2895e76c431SBarry Smith { 290393d2d9aSLois Curfman McInnes int ierr; 2913a40ed3dSBarry Smith 2923a40ed3dSBarry Smith PetscFunctionBegin; 2935baf8537SBarry Smith if (snes->nwork) { 2944b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 29521c89e3eSBarry Smith } 2965c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 2985e76c431SBarry Smith } 29904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3004a2ae208SSatish Balay #undef __FUNCT__ 3014a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch" 30282bf6240SBarry Smith 3034b828684SBarry Smith /*@C 3045e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 3055e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3065e76c431SBarry Smith to serve as a template and is not recommended for general use. 3075e76c431SBarry Smith 308c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 309c7afd0dbSLois Curfman McInnes 3105e76c431SBarry Smith Input Parameters: 311c7afd0dbSLois Curfman McInnes + snes - nonlinear context 312af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3135e76c431SBarry Smith . x - current iterate 3145e76c431SBarry Smith . f - residual evaluated at x 3155e76c431SBarry Smith . y - search direction (contains new iterate on output) 3165e76c431SBarry Smith . w - work vector 317c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3185e76c431SBarry Smith 319c4a48953SLois Curfman McInnes Output Parameters: 320c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3215e76c431SBarry Smith . y - new iterate (contains search direction on input) 3225e76c431SBarry Smith . gnorm - 2-norm of g 3235e76c431SBarry Smith . ynorm - 2-norm of search length 324c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 325fee21e36SBarry Smith 326c4a48953SLois Curfman McInnes Options Database Key: 3274b27c08aSLois Curfman McInnes . -snes_ls basic - Activates SNESNoLineSearch() 328c4a48953SLois Curfman McInnes 32936851e7fSLois Curfman McInnes Level: advanced 33036851e7fSLois Curfman McInnes 33128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 33228ae5a21SLois Curfman McInnes 333f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 334af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3355e76c431SBarry Smith @*/ 3365d1a10b1SBarry 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) 3375e76c431SBarry Smith { 3385e42470aSBarry Smith int ierr; 339ea709b57SSatish Balay PetscScalar mone = -1.0; 3404b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3418f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3425334005bSBarry Smith 3433a40ed3dSBarry Smith PetscFunctionBegin; 344761aaf1bSLois Curfman McInnes *flag = 0; 345d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 346a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 347ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3488f99978dSLois Curfman McInnes if (neP->CheckStep) { 3498f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3508f99978dSLois Curfman McInnes } 351ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 352a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 353d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3543a40ed3dSBarry Smith PetscFunctionReturn(0); 3555e76c431SBarry Smith } 35604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3574a2ae208SSatish Balay #undef __FUNCT__ 3584a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms" 35982bf6240SBarry Smith 36029e0b56fSBarry Smith /*@C 36129e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 36229e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 36329e0b56fSBarry Smith even compute the norm of the function or search direction; this 36429e0b56fSBarry Smith is intended only when you know the full step is fine and are 36529e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 366c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 367c7afd0dbSLois Curfman McInnes 368c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 36929e0b56fSBarry Smith 37029e0b56fSBarry Smith Input Parameters: 371c7afd0dbSLois Curfman McInnes + snes - nonlinear context 372af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 37329e0b56fSBarry Smith . x - current iterate 37429e0b56fSBarry Smith . f - residual evaluated at x 37529e0b56fSBarry Smith . y - search direction (contains new iterate on output) 37629e0b56fSBarry Smith . w - work vector 377c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 37829e0b56fSBarry Smith 37929e0b56fSBarry Smith Output Parameters: 380c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 38129e0b56fSBarry Smith . gnorm - not changed 38229e0b56fSBarry Smith . ynorm - not changed 383c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 384fee21e36SBarry Smith 38529e0b56fSBarry Smith Options Database Key: 3864b27c08aSLois Curfman McInnes . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 38729e0b56fSBarry Smith 3888cbba510SBarry Smith Notes: 389ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 390ea56f5baSLois Curfman McInnes either the options 391ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 392ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 393329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 394329f5518SBarry Smith otherwise, the SNES solver will generate an error. 395329f5518SBarry Smith 396329f5518SBarry Smith During the final iteration this will not evaluate the function at 397329f5518SBarry Smith the solution point. This is to save a function evaluation while 398329f5518SBarry Smith using pseudo-timestepping. 3998cbba510SBarry Smith 400ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 401ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 402ea56f5baSLois Curfman McInnes correct, since they are not computed. 403ea56f5baSLois Curfman McInnes 404ea56f5baSLois Curfman McInnes Level: advanced 4058cbba510SBarry Smith 40629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 40729e0b56fSBarry Smith 40829e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 40929e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 41029e0b56fSBarry Smith @*/ 4115d1a10b1SBarry 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) 41229e0b56fSBarry Smith { 41329e0b56fSBarry Smith int ierr; 414ea709b57SSatish Balay PetscScalar mone = -1.0; 4154b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4168f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 41729e0b56fSBarry Smith 4183a40ed3dSBarry Smith PetscFunctionBegin; 4198cbba510SBarry Smith *flag = 0; 420d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 42129e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4228f99978dSLois Curfman McInnes if (neP->CheckStep) { 4238f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4248f99978dSLois Curfman McInnes } 425329f5518SBarry Smith 426329f5518SBarry Smith /* don't evaluate function the last time through */ 427329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 42829e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 429329f5518SBarry Smith } 430d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4313a40ed3dSBarry Smith PetscFunctionReturn(0); 43229e0b56fSBarry Smith } 43304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4344a2ae208SSatish Balay #undef __FUNCT__ 4354a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch" 4364b828684SBarry Smith /*@C 437f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4385e76c431SBarry Smith 439c7afd0dbSLois Curfman McInnes Collective on SNES 440c7afd0dbSLois Curfman McInnes 4415e76c431SBarry Smith Input Parameters: 442c7afd0dbSLois Curfman McInnes + snes - nonlinear context 443af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4445e76c431SBarry Smith . x - current iterate 4455e76c431SBarry Smith . f - residual evaluated at x 4465e76c431SBarry Smith . y - search direction (contains new iterate on output) 4475e76c431SBarry Smith . w - work vector 448c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4495e76c431SBarry Smith 450393d2d9aSLois Curfman McInnes Output Parameters: 451c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4525e76c431SBarry Smith . y - new iterate (contains search direction on input) 4535e76c431SBarry Smith . gnorm - 2-norm of g 4545e76c431SBarry Smith . ynorm - 2-norm of search length 455c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 456fee21e36SBarry Smith 457c4a48953SLois Curfman McInnes Options Database Key: 4584b27c08aSLois Curfman McInnes $ -snes_ls cubic - Activates SNESCubicLineSearch() 459c4a48953SLois Curfman McInnes 4605e76c431SBarry Smith Notes: 4615e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4625e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4635e76c431SBarry Smith 46436851e7fSLois Curfman McInnes Level: advanced 46536851e7fSLois Curfman McInnes 46628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 46728ae5a21SLois Curfman McInnes 468af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4695e76c431SBarry Smith @*/ 4705d1a10b1SBarry 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) 4715e76c431SBarry Smith { 472406556e6SLois Curfman McInnes /* 473406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 474406556e6SLois Curfman McInnes minimization problem: 475406556e6SLois Curfman McInnes min z(x): R^n -> R, 476406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 477406556e6SLois Curfman McInnes */ 478406556e6SLois Curfman McInnes 4795ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 480329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 481aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 482ea709b57SSatish Balay PetscScalar cinitslope,clambda; 4836b5873e3SBarry Smith #endif 4845e42470aSBarry Smith int ierr,count; 4854b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 486ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 4878f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4885e76c431SBarry Smith 4893a40ed3dSBarry Smith PetscFunctionBegin; 490d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 491761aaf1bSLois Curfman McInnes *flag = 0; 4925e76c431SBarry Smith alpha = neP->alpha; 4935e76c431SBarry Smith maxstep = neP->maxstep; 4945e76c431SBarry Smith steptol = neP->steptol; 4955e76c431SBarry Smith 496cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 49763b9fa5eSBarry Smith if (*ynorm == 0.0) { 49863b9fa5eSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n"); 499a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 500a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 501a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 50263b9fa5eSBarry Smith *flag = -1; 503ad922ac9SBarry Smith goto theend1; 50494a9d846SBarry Smith } 5055e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5065e42470aSBarry Smith scale = maxstep/(*ynorm); 507aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 508b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 5096b5873e3SBarry Smith #else 510b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 5116b5873e3SBarry Smith #endif 512393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5135e76c431SBarry Smith *ynorm = maxstep; 5145e76c431SBarry Smith } 515ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5165ca10a37SBarry Smith minlambda = steptol/rellength; 517a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 518aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 519a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 520329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5215e42470aSBarry Smith #else 522a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5235e42470aSBarry Smith #endif 5245e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5255e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5265e76c431SBarry Smith 527393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5285334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 52978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 530313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5315d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 532393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 533b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 534ad922ac9SBarry Smith goto theend1; 5355e76c431SBarry Smith } 5365e76c431SBarry Smith 5375e76c431SBarry Smith /* Fit points with quadratic */ 538313b5538SBarry Smith lambda = 1.0; 539a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5405e76c431SBarry Smith lambdaprev = lambda; 5415e76c431SBarry Smith gnormprev = *gnorm; 542329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 543ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 544ddd12767SBarry Smith else lambda = lambdatemp; 545393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 546ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 547aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 548ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5495e42470aSBarry Smith #else 550ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5515e42470aSBarry Smith #endif 55278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 553cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5545ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 555393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 5565ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda); 557ad922ac9SBarry Smith goto theend1; 5585e76c431SBarry Smith } 5595e76c431SBarry Smith 5605e76c431SBarry Smith /* Fit points with cubic */ 5615e76c431SBarry Smith count = 1; 5628229bfc2SMatthew Knepley while (count < 10000) { 5635e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 564b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 5655ed2d596SBarry 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); 566f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 567761aaf1bSLois Curfman McInnes *flag = -1; break; 5685e76c431SBarry Smith } 5695d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 5705d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 571ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5722b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5735e76c431SBarry Smith d = b*b - 3*a*initslope; 5745e76c431SBarry Smith if (d < 0.0) d = 0.0; 5755e76c431SBarry Smith if (a == 0.0) { 5765e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5775e76c431SBarry Smith } else { 5785e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5795e76c431SBarry Smith } 5805e76c431SBarry Smith lambdaprev = lambda; 5815e76c431SBarry Smith gnormprev = *gnorm; 582329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 583329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5845e76c431SBarry Smith else lambda = lambdatemp; 585393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 586ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 587aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 588ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 589393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5905e42470aSBarry Smith #else 591ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5925e42470aSBarry Smith #endif 59378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 594cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5955ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 596393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 5975ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda); 5982715565aSLois Curfman McInnes goto theend1; 5992b022350SLois Curfman McInnes } else { 6005ed2d596SBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda); 6015e76c431SBarry Smith } 6025e76c431SBarry Smith count++; 6035e76c431SBarry Smith } 6048229bfc2SMatthew Knepley if (count >= 10000) { 605cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6068229bfc2SMatthew Knepley } 607ad922ac9SBarry Smith theend1: 6088f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6098f99978dSLois Curfman McInnes if (neP->CheckStep) { 6108f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 6118f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6128f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6138f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6148f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6158f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6168f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6178f99978dSLois Curfman McInnes } 6188f99978dSLois Curfman McInnes } 619d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6203a40ed3dSBarry Smith PetscFunctionReturn(0); 6215e76c431SBarry Smith } 62204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6234a2ae208SSatish Balay #undef __FUNCT__ 6244a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch" 6254b828684SBarry Smith /*@C 626f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6275e76c431SBarry Smith 628c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 629c7afd0dbSLois Curfman McInnes 6305e42470aSBarry Smith Input Parameters: 631c7afd0dbSLois Curfman McInnes + snes - the SNES context 632af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6335e42470aSBarry Smith . x - current iterate 6345e42470aSBarry Smith . f - residual evaluated at x 6355e42470aSBarry Smith . y - search direction (contains new iterate on output) 6365e42470aSBarry Smith . w - work vector 637c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6385e42470aSBarry Smith 639c4a48953SLois Curfman McInnes Output Parameters: 640c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6415e42470aSBarry Smith . y - new iterate (contains search direction on input) 6425e42470aSBarry Smith . gnorm - 2-norm of g 6435e42470aSBarry Smith . ynorm - 2-norm of search length 644c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 645fee21e36SBarry Smith 646c4a48953SLois Curfman McInnes Options Database Key: 6474b27c08aSLois Curfman McInnes . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 648c4a48953SLois Curfman McInnes 6495e42470aSBarry Smith Notes: 6504b27c08aSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNESLS method. 6515e42470aSBarry Smith 65236851e7fSLois Curfman McInnes Level: advanced 65336851e7fSLois Curfman McInnes 654f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 655f59ffdedSLois Curfman McInnes 656af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6575e42470aSBarry Smith @*/ 6585d1a10b1SBarry 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) 6595e76c431SBarry Smith { 660406556e6SLois Curfman McInnes /* 661406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 662406556e6SLois Curfman McInnes minimization problem: 663406556e6SLois Curfman McInnes min z(x): R^n -> R, 664406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 665406556e6SLois Curfman McInnes */ 6665ca10a37SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 667aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 668ea709b57SSatish Balay PetscScalar cinitslope,clambda; 6696b5873e3SBarry Smith #endif 6705e42470aSBarry Smith int ierr,count; 6714b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 672ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 6738f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6745e76c431SBarry Smith 6753a40ed3dSBarry Smith PetscFunctionBegin; 676d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 677761aaf1bSLois Curfman McInnes *flag = 0; 6785e42470aSBarry Smith alpha = neP->alpha; 6795e42470aSBarry Smith maxstep = neP->maxstep; 6805e42470aSBarry Smith steptol = neP->steptol; 6815e76c431SBarry Smith 6823505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 68363b9fa5eSBarry Smith if (*ynorm == 0.0) { 684b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 685b37302c6SLois Curfman McInnes *gnorm = fnorm; 686b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 687b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 68863b9fa5eSBarry Smith *flag = -1; 689ad922ac9SBarry Smith goto theend2; 69094a9d846SBarry Smith } 6915e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6925e42470aSBarry Smith scale = maxstep/(*ynorm); 693393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6945e42470aSBarry Smith *ynorm = maxstep; 6955e76c431SBarry Smith } 696ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 6975ca10a37SBarry Smith minlambda = steptol/rellength; 698a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 699aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 700a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 701329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7025e42470aSBarry Smith #else 703a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7045e42470aSBarry Smith #endif 7055e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7065e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7075e42470aSBarry Smith 708393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 7095334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 71078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 711cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7125d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 713393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 714b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 715ad922ac9SBarry Smith goto theend2; 7165e42470aSBarry Smith } 7175e42470aSBarry Smith 7185e42470aSBarry Smith /* Fit points with quadratic */ 719313b5538SBarry Smith lambda = 1.0; 7205e42470aSBarry Smith count = 1; 7215ca10a37SBarry Smith while (PETSC_TRUE) { 7225e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 723b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 724b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 725f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 726761aaf1bSLois Curfman McInnes *flag = -1; break; 7275e42470aSBarry Smith } 728a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 729329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 730329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 731329f5518SBarry Smith else lambda = lambdatemp; 732393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7333505fcc1SBarry Smith lambdaneg = -lambda; 734aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7353505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7365e42470aSBarry Smith #else 7373505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7385e42470aSBarry Smith #endif 73978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 740cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7415ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 742393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 743b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 74406259719SBarry Smith break; 7455e42470aSBarry Smith } 7465e42470aSBarry Smith count++; 7475e42470aSBarry Smith } 748ad922ac9SBarry Smith theend2: 7498f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7508f99978dSLois Curfman McInnes if (neP->CheckStep) { 7518f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7528f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7538f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7548f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7558f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7568f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7578f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7588f99978dSLois Curfman McInnes } 7598f99978dSLois Curfman McInnes } 760d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7613a40ed3dSBarry Smith PetscFunctionReturn(0); 7625e76c431SBarry Smith } 76304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7644a2ae208SSatish Balay #undef __FUNCT__ 7654a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 766c9e6a524SLois Curfman McInnes /*@C 76777c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7684b27c08aSLois Curfman McInnes by the method SNESLS. 7695e76c431SBarry Smith 7705e76c431SBarry Smith Input Parameters: 771c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 772af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 773c7afd0dbSLois Curfman McInnes - func - pointer to int function 7745e76c431SBarry Smith 775fee21e36SBarry Smith Collective on SNES 776fee21e36SBarry Smith 777c4a48953SLois Curfman McInnes Available Routines: 778c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 779f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 780f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 781af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7825e76c431SBarry Smith 783c4a48953SLois Curfman McInnes Options Database Keys: 7844b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 7854b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 7864b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 7874b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 7883304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 7893304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 790c4a48953SLois Curfman McInnes 7915e76c431SBarry Smith Calling sequence of func: 792bd208895SLois Curfman McInnes .vb 793af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 794329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 795bd208895SLois Curfman McInnes .ve 7965e76c431SBarry Smith 7975e76c431SBarry Smith Input parameters for func: 798c7afd0dbSLois Curfman McInnes + snes - nonlinear context 799af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8005e76c431SBarry Smith . x - current iterate 8015e76c431SBarry Smith . f - residual evaluated at x 8025e76c431SBarry Smith . y - search direction (contains new iterate on output) 8035e76c431SBarry Smith . w - work vector 804c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 8055e76c431SBarry Smith 8065e76c431SBarry Smith Output parameters for func: 807c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 8085e76c431SBarry Smith . y - new iterate (contains search direction on input) 8095e76c431SBarry Smith . gnorm - 2-norm of g 8105e76c431SBarry Smith . ynorm - 2-norm of search length 811c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 812761aaf1bSLois Curfman McInnes on failure. 813f59ffdedSLois Curfman McInnes 81436851e7fSLois Curfman McInnes Level: advanced 81536851e7fSLois Curfman McInnes 816f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 817f59ffdedSLois Curfman McInnes 8185d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 8195d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 8205e76c431SBarry Smith @*/ 821329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 8225e76c431SBarry Smith { 823329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 82482bf6240SBarry Smith 8253a40ed3dSBarry Smith PetscFunctionBegin; 826c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 82782bf6240SBarry Smith if (f) { 828af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 82982bf6240SBarry Smith } 8303a40ed3dSBarry Smith PetscFunctionReturn(0); 8315e76c431SBarry Smith } 83204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 833fb2e594dSBarry Smith EXTERN_C_BEGIN 8344a2ae208SSatish Balay #undef __FUNCT__ 8354a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS" 836af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 83787828ca2SBarry Smith PetscReal,PetscReal*,PetscReal*,int*),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 /* -------------------------------------------------------------------------- */ 910c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 9114a2ae208SSatish Balay #undef __FUNCT__ 9124a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS" 9138f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 914c8dd0c0dSLois Curfman McInnes { 915c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 9164b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 9174b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 918c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 919c8dd0c0dSLois Curfman McInnes } 920c8dd0c0dSLois Curfman McInnes EXTERN_C_END 921c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92204d965bbSLois Curfman McInnes /* 9234b27c08aSLois Curfman McInnes SNESPrintHelp_LS - Prints all options for the SNES_LS method. 924329e7e40SMatthew Knepley 925329e7e40SMatthew Knepley Input Parameter: 926329e7e40SMatthew Knepley . snes - the SNES context 927329e7e40SMatthew Knepley 928329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 929329e7e40SMatthew Knepley */ 930329e7e40SMatthew Knepley #undef __FUNCT__ 9314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS" 9324b27c08aSLois Curfman McInnes static int SNESPrintHelp_LS(SNES snes,char *p) 933329e7e40SMatthew Knepley { 9344b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 935329e7e40SMatthew Knepley 936329e7e40SMatthew Knepley PetscFunctionBegin; 9374b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 9384b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 9394b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 9404b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 9414b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 942329e7e40SMatthew Knepley PetscFunctionReturn(0); 943329e7e40SMatthew Knepley } 944329e7e40SMatthew Knepley 945329e7e40SMatthew Knepley /* 9464b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 94704d965bbSLois Curfman McInnes 94804d965bbSLois Curfman McInnes Input Parameters: 94904d965bbSLois Curfman McInnes . SNES - the SNES context 95004d965bbSLois Curfman McInnes . viewer - visualization context 95104d965bbSLois Curfman McInnes 95204d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 95304d965bbSLois Curfman McInnes */ 9544a2ae208SSatish Balay #undef __FUNCT__ 9554b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 9564b27c08aSLois Curfman McInnes static int SNESView_LS(SNES snes,PetscViewer viewer) 957a935fc98SLois Curfman McInnes { 9584b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 95919bcc07fSBarry Smith char *cstr; 96051695f54SLois Curfman McInnes int ierr; 9616831982aSBarry Smith PetscTruth isascii; 962a935fc98SLois Curfman McInnes 9633a40ed3dSBarry Smith PetscFunctionBegin; 964b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 9650f5bd95cSBarry Smith if (isascii) { 96619bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 96719bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 96819bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 96919bcc07fSBarry Smith else cstr = "unknown"; 970b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 971b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9725cd90555SBarry Smith } else { 97329bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 97419bcc07fSBarry Smith } 9753a40ed3dSBarry Smith PetscFunctionReturn(0); 976a935fc98SLois Curfman McInnes } 97704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 97804d965bbSLois Curfman McInnes /* 9794b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 98004d965bbSLois Curfman McInnes 98104d965bbSLois Curfman McInnes Input Parameter: 98204d965bbSLois Curfman McInnes . snes - the SNES context 98304d965bbSLois Curfman McInnes 98404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 98504d965bbSLois Curfman McInnes */ 9864a2ae208SSatish Balay #undef __FUNCT__ 9874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 9884b27c08aSLois Curfman McInnes static int SNESSetFromOptions_LS(SNES snes) 9895e42470aSBarry Smith { 9904b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 991b49fd9e1SBarry Smith char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 99222e36657SBarry Smith int ierr,indx; 993f1af5d2fSBarry Smith PetscTruth flg; 9945e42470aSBarry Smith 9953a40ed3dSBarry Smith PetscFunctionBegin; 996b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 9974b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 9984b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 9994b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1000186905e3SBarry Smith 100122e36657SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 100225cdf11fSBarry Smith if (flg) { 100322e36657SBarry Smith switch (indx) { 1004b49fd9e1SBarry Smith case 0: 1005af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1006b49fd9e1SBarry Smith break; 1007b49fd9e1SBarry Smith case 1: 1008af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1009b49fd9e1SBarry Smith break; 1010b49fd9e1SBarry Smith case 2: 1011af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1012b49fd9e1SBarry Smith break; 1013b49fd9e1SBarry Smith case 3: 1014af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1015b49fd9e1SBarry Smith break; 10165e42470aSBarry Smith } 10175e42470aSBarry Smith } 1018b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10193a40ed3dSBarry Smith PetscFunctionReturn(0); 10205e42470aSBarry Smith } 102104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 102204d965bbSLois Curfman McInnes /* 10234b27c08aSLois Curfman McInnes SNESCreate_LS - Creates a nonlinear solver context for the SNESLS method, 10244b27c08aSLois Curfman McInnes SNES_LS, and sets this as the private data within the generic nonlinear solver 102504d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 102604d965bbSLois Curfman McInnes 102704d965bbSLois Curfman McInnes 102804d965bbSLois Curfman McInnes Input Parameter: 102904d965bbSLois Curfman McInnes . snes - the SNES context 103004d965bbSLois Curfman McInnes 103104d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 103204d965bbSLois Curfman McInnes */ 1033fb2e594dSBarry Smith EXTERN_C_BEGIN 10344a2ae208SSatish Balay #undef __FUNCT__ 10354b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 10364b27c08aSLois Curfman McInnes int SNESCreate_LS(SNES snes) 10375e42470aSBarry Smith { 103882bf6240SBarry Smith int ierr; 10394b27c08aSLois Curfman McInnes SNES_LS *neP; 10405e42470aSBarry Smith 10413a40ed3dSBarry Smith PetscFunctionBegin; 10424b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 10434b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 10444b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 10454b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 10464b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 10474b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 10484b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 10495baf8537SBarry Smith snes->nwork = 0; 10505e42470aSBarry Smith 10514b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 10524b27c08aSLois Curfman McInnes PetscLogObjectMemory(snes,sizeof(SNES_LS)); 10535e42470aSBarry Smith snes->data = (void*)neP; 10545e42470aSBarry Smith neP->alpha = 1.e-4; 10555e42470aSBarry Smith neP->maxstep = 1.e8; 10565e42470aSBarry Smith neP->steptol = 1.e-12; 10575e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1058c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1059c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1060c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 106182bf6240SBarry Smith 10625d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10635d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 106482bf6240SBarry Smith 10653a40ed3dSBarry Smith PetscFunctionReturn(0); 10665e42470aSBarry Smith } 1067fb2e594dSBarry Smith EXTERN_C_END 106882bf6240SBarry Smith 106982bf6240SBarry Smith 107082bf6240SBarry Smith 1071