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); 26b0a32e0cSBarry Smith PetscLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 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); 40b0a32e0cSBarry Smith PetscLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 4774637425SBarry Smith Checks if J^T(F - AX) = 0 4874637425SBarry Smith */ 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) { 69b0a32e0cSBarry Smith PetscLogInfo(0,"SNESSolve_EQ_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 8804d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 8904d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 9004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9304d965bbSLois Curfman McInnes for all nonlinear solvers. 9404d965bbSLois Curfman McInnes 9504d965bbSLois Curfman McInnes Another key routine is: 9604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9904d965bbSLois Curfman McInnes SNESSolve() if necessary. 10004d965bbSLois Curfman McInnes 10104d965bbSLois Curfman McInnes Additional basic routines are: 10204d965bbSLois Curfman McInnes 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 /* 11304d965bbSLois Curfman McInnes SNESSolve_EQ_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__ 1324a2ae208SSatish Balay #define __FUNCT__ "SNESSolve_EQ_LS" 133f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 1345e76c431SBarry Smith { 1356831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_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; 1405e76c431SBarry Smith 1413a40ed3dSBarry Smith PetscFunctionBegin; 142184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 143184914b5SBarry Smith 1445e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1455e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14639e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1475e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1485e42470aSBarry Smith G = snes->work[1]; 1495e42470aSBarry Smith W = snes->work[2]; 1505e76c431SBarry Smith 1514c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1524c49b128SBarry Smith snes->iter = 0; 1534c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1545334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 155cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1560f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1575e42470aSBarry Smith snes->norm = fnorm; 1580f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15984c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16094a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1615e76c431SBarry Smith 162184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16394a9d846SBarry Smith 164d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 165d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 166d034289bSBarry Smith 1675e76c431SBarry Smith for (i=0; i<maxits; i++) { 1685e76c431SBarry Smith 169*329e7e40SMatthew Knepley /* Call general purpose update function */ 170*329e7e40SMatthew Knepley ierr = (*snes->update)(snes, snes->iter); CHKERRQ(ierr); 171*329e7e40SMatthew Knepley 172ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1735334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1745334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 17578b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 17674637425SBarry Smith 177b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 17874637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 17985471664SBarry Smith } 18074637425SBarry Smith 18190cbd584SBarry Smith /* should check what happened to the linear solve? */ 1823505fcc1SBarry Smith snes->linear_its += lits; 183b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 184ea4d3ed3SLois Curfman McInnes 185ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 186ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 187ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 188ea4d3ed3SLois Curfman McInnes */ 18981b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 190af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 191b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 192a3b891d8SBarry Smith 193a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 194a3b891d8SBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 195a3b891d8SBarry Smith fnorm = gnorm; 196a3b891d8SBarry Smith 1973505fcc1SBarry Smith if (lsfail) { 1988a5d9ee4SBarry Smith PetscTruth ismin; 1993505fcc1SBarry Smith snes->nfailures++; 2003505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2018a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2028a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2033505fcc1SBarry Smith break; 2043505fcc1SBarry Smith } 2055e76c431SBarry Smith 2060f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 20784c569c9SBarry Smith snes->iter = i+1; 2085e42470aSBarry Smith snes->norm = fnorm; 2090f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 21084c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 21194a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 2125e76c431SBarry Smith 2135e76c431SBarry Smith /* Test for convergence */ 21429e0b56fSBarry Smith if (snes->converged) { 21529e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2163505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2173505fcc1SBarry Smith if (snes->reason) { 21816c95cacSBarry Smith break; 21916c95cacSBarry Smith } 22016c95cacSBarry Smith } 22129e0b56fSBarry Smith } 22239e2f89bSBarry Smith if (X != snes->vec_sol) { 223393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 224e15f7bb5SBarry Smith } 225e15f7bb5SBarry Smith if (F != snes->vec_func) { 226e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 227e15f7bb5SBarry Smith } 22839e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 22939e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 23052392280SLois Curfman McInnes if (i == maxits) { 231b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 23252392280SLois Curfman McInnes i--; 2333505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 23452392280SLois Curfman McInnes } 2350f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2360f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2375e42470aSBarry Smith *outits = i+1; 2383a40ed3dSBarry Smith PetscFunctionReturn(0); 2395e76c431SBarry Smith } 24004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 24104d965bbSLois Curfman McInnes /* 24204d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 2436831982aSBarry Smith of the SNESEQLS nonlinear solver. 24404d965bbSLois Curfman McInnes 24504d965bbSLois Curfman McInnes Input Parameter: 24604d965bbSLois Curfman McInnes . snes - the SNES context 24704d965bbSLois Curfman McInnes . x - the solution vector 24804d965bbSLois Curfman McInnes 24904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 25004d965bbSLois Curfman McInnes 25104d965bbSLois Curfman McInnes Notes: 25204d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 25304d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 25404d965bbSLois Curfman McInnes the call to SNESSolve(). 25504d965bbSLois Curfman McInnes */ 2564a2ae208SSatish Balay #undef __FUNCT__ 2574a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_LS" 258f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 2595e76c431SBarry Smith { 2605e42470aSBarry Smith int ierr; 2613a40ed3dSBarry Smith 2623a40ed3dSBarry Smith PetscFunctionBegin; 26381b6cf68SLois Curfman McInnes snes->nwork = 4; 264d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 265b0a32e0cSBarry Smith PetscLogObjectParents(snes,snes->nwork,snes->work); 26681b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2685e76c431SBarry Smith } 26904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27004d965bbSLois Curfman McInnes /* 2716831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 27204d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 27304d965bbSLois Curfman McInnes 27404d965bbSLois Curfman McInnes Input Parameter: 27504d965bbSLois Curfman McInnes . snes - the SNES context 27604d965bbSLois Curfman McInnes 277de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 27804d965bbSLois Curfman McInnes */ 2794a2ae208SSatish Balay #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_LS" 281e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2825e76c431SBarry Smith { 283393d2d9aSLois Curfman McInnes int ierr; 2843a40ed3dSBarry Smith 2853a40ed3dSBarry Smith PetscFunctionBegin; 2865baf8537SBarry Smith if (snes->nwork) { 2874b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 28821c89e3eSBarry Smith } 2895c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2903a40ed3dSBarry Smith PetscFunctionReturn(0); 2915e76c431SBarry Smith } 29204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2934a2ae208SSatish Balay #undef __FUNCT__ 2944a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch" 29582bf6240SBarry Smith 2964b828684SBarry Smith /*@C 2975e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2985e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2995e76c431SBarry Smith to serve as a template and is not recommended for general use. 3005e76c431SBarry Smith 301c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 302c7afd0dbSLois Curfman McInnes 3035e76c431SBarry Smith Input Parameters: 304c7afd0dbSLois Curfman McInnes + snes - nonlinear context 305af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3065e76c431SBarry Smith . x - current iterate 3075e76c431SBarry Smith . f - residual evaluated at x 3085e76c431SBarry Smith . y - search direction (contains new iterate on output) 3095e76c431SBarry Smith . w - work vector 310c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3115e76c431SBarry Smith 312c4a48953SLois Curfman McInnes Output Parameters: 313c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3145e76c431SBarry Smith . y - new iterate (contains search direction on input) 3155e76c431SBarry Smith . gnorm - 2-norm of g 3165e76c431SBarry Smith . ynorm - 2-norm of search length 317c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 318fee21e36SBarry Smith 319c4a48953SLois Curfman McInnes Options Database Key: 320c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 321c4a48953SLois Curfman McInnes 32236851e7fSLois Curfman McInnes Level: advanced 32336851e7fSLois Curfman McInnes 32428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 32528ae5a21SLois Curfman McInnes 326f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 327af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3285e76c431SBarry Smith @*/ 3295d1a10b1SBarry 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) 3305e76c431SBarry Smith { 3315e42470aSBarry Smith int ierr; 332ea709b57SSatish Balay PetscScalar mone = -1.0; 3336831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3348f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3355334005bSBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 337761aaf1bSLois Curfman McInnes *flag = 0; 338b0a32e0cSBarry Smith ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 339a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 340ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3418f99978dSLois Curfman McInnes if (neP->CheckStep) { 3428f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3438f99978dSLois Curfman McInnes } 344ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 345a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 346b0a32e0cSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 3485e76c431SBarry Smith } 34904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3504a2ae208SSatish Balay #undef __FUNCT__ 3514a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms" 35282bf6240SBarry Smith 35329e0b56fSBarry Smith /*@C 35429e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 35529e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 35629e0b56fSBarry Smith even compute the norm of the function or search direction; this 35729e0b56fSBarry Smith is intended only when you know the full step is fine and are 35829e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 359c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 360c7afd0dbSLois Curfman McInnes 361c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 36229e0b56fSBarry Smith 36329e0b56fSBarry Smith Input Parameters: 364c7afd0dbSLois Curfman McInnes + snes - nonlinear context 365af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 36629e0b56fSBarry Smith . x - current iterate 36729e0b56fSBarry Smith . f - residual evaluated at x 36829e0b56fSBarry Smith . y - search direction (contains new iterate on output) 36929e0b56fSBarry Smith . w - work vector 370c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 37129e0b56fSBarry Smith 37229e0b56fSBarry Smith Output Parameters: 373c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 37429e0b56fSBarry Smith . gnorm - not changed 37529e0b56fSBarry Smith . ynorm - not changed 376c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 377fee21e36SBarry Smith 37829e0b56fSBarry Smith Options Database Key: 379c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 38029e0b56fSBarry Smith 3818cbba510SBarry Smith Notes: 382ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 383ea56f5baSLois Curfman McInnes either the options 384ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 385ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 386329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 387329f5518SBarry Smith otherwise, the SNES solver will generate an error. 388329f5518SBarry Smith 389329f5518SBarry Smith During the final iteration this will not evaluate the function at 390329f5518SBarry Smith the solution point. This is to save a function evaluation while 391329f5518SBarry Smith using pseudo-timestepping. 3928cbba510SBarry Smith 393ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 394ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 395ea56f5baSLois Curfman McInnes correct, since they are not computed. 396ea56f5baSLois Curfman McInnes 397ea56f5baSLois Curfman McInnes Level: advanced 3988cbba510SBarry Smith 39929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 40029e0b56fSBarry Smith 40129e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 40229e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 40329e0b56fSBarry Smith @*/ 4045d1a10b1SBarry 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) 40529e0b56fSBarry Smith { 40629e0b56fSBarry Smith int ierr; 407ea709b57SSatish Balay PetscScalar mone = -1.0; 4086831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4098f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 41029e0b56fSBarry Smith 4113a40ed3dSBarry Smith PetscFunctionBegin; 4128cbba510SBarry Smith *flag = 0; 413b0a32e0cSBarry Smith ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 41429e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4158f99978dSLois Curfman McInnes if (neP->CheckStep) { 4168f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4178f99978dSLois Curfman McInnes } 418329f5518SBarry Smith 419329f5518SBarry Smith /* don't evaluate function the last time through */ 420329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 42129e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 422329f5518SBarry Smith } 423b0a32e0cSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4243a40ed3dSBarry Smith PetscFunctionReturn(0); 42529e0b56fSBarry Smith } 42604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4274a2ae208SSatish Balay #undef __FUNCT__ 4284a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch" 4294b828684SBarry Smith /*@C 430f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4315e76c431SBarry Smith 432c7afd0dbSLois Curfman McInnes Collective on SNES 433c7afd0dbSLois Curfman McInnes 4345e76c431SBarry Smith Input Parameters: 435c7afd0dbSLois Curfman McInnes + snes - nonlinear context 436af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4375e76c431SBarry Smith . x - current iterate 4385e76c431SBarry Smith . f - residual evaluated at x 4395e76c431SBarry Smith . y - search direction (contains new iterate on output) 4405e76c431SBarry Smith . w - work vector 441c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4425e76c431SBarry Smith 443393d2d9aSLois Curfman McInnes Output Parameters: 444c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4455e76c431SBarry Smith . y - new iterate (contains search direction on input) 4465e76c431SBarry Smith . gnorm - 2-norm of g 4475e76c431SBarry Smith . ynorm - 2-norm of search length 448c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 449fee21e36SBarry Smith 450c4a48953SLois Curfman McInnes Options Database Key: 451c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 452c4a48953SLois Curfman McInnes 4535e76c431SBarry Smith Notes: 4545e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4555e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4565e76c431SBarry Smith 45736851e7fSLois Curfman McInnes Level: advanced 45836851e7fSLois Curfman McInnes 45928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 46028ae5a21SLois Curfman McInnes 461af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4625e76c431SBarry Smith @*/ 4635d1a10b1SBarry 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) 4645e76c431SBarry Smith { 465406556e6SLois Curfman McInnes /* 466406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 467406556e6SLois Curfman McInnes minimization problem: 468406556e6SLois Curfman McInnes min z(x): R^n -> R, 469406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 470406556e6SLois Curfman McInnes */ 471406556e6SLois Curfman McInnes 472329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 473329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 474aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 475ea709b57SSatish Balay PetscScalar cinitslope,clambda; 4766b5873e3SBarry Smith #endif 4775e42470aSBarry Smith int ierr,count; 4786831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 479ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 4808f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4815e76c431SBarry Smith 4823a40ed3dSBarry Smith PetscFunctionBegin; 483b0a32e0cSBarry Smith ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 484761aaf1bSLois Curfman McInnes *flag = 0; 4855e76c431SBarry Smith alpha = neP->alpha; 4865e76c431SBarry Smith maxstep = neP->maxstep; 4875e76c431SBarry Smith steptol = neP->steptol; 4885e76c431SBarry Smith 489cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 490a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 491b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 492a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 493a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 494a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 495ad922ac9SBarry Smith goto theend1; 49694a9d846SBarry Smith } 4975e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4985e42470aSBarry Smith scale = maxstep/(*ynorm); 499aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 500b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 5016b5873e3SBarry Smith #else 502b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 5036b5873e3SBarry Smith #endif 504393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5055e76c431SBarry Smith *ynorm = maxstep; 5065e76c431SBarry Smith } 5075e76c431SBarry Smith minlambda = steptol/(*ynorm); 508a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 509aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 510a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 511329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5125e42470aSBarry Smith #else 513a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5145e42470aSBarry Smith #endif 5155e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5165e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5175e76c431SBarry Smith 518393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5195334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 52078b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 521313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5225d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 523393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 524b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 525ad922ac9SBarry Smith goto theend1; 5265e76c431SBarry Smith } 5275e76c431SBarry Smith 5285e76c431SBarry Smith /* Fit points with quadratic */ 529313b5538SBarry Smith lambda = 1.0; 530a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5315e76c431SBarry Smith lambdaprev = lambda; 5325e76c431SBarry Smith gnormprev = *gnorm; 533329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 534ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 535ddd12767SBarry Smith else lambda = lambdatemp; 536393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 537ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 538aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 539ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5405e42470aSBarry Smith #else 541ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5425e42470aSBarry Smith #endif 54378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 544cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5455d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 546393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 547b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 548ad922ac9SBarry Smith goto theend1; 5495e76c431SBarry Smith } 5505e76c431SBarry Smith 5515e76c431SBarry Smith /* Fit points with cubic */ 5525e76c431SBarry Smith count = 1; 5535e76c431SBarry Smith while (1) { 5545e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 555b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 556b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 557393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 558761aaf1bSLois Curfman McInnes *flag = -1; break; 5595e76c431SBarry Smith } 5605d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 5615d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 562ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5632b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5645e76c431SBarry Smith d = b*b - 3*a*initslope; 5655e76c431SBarry Smith if (d < 0.0) d = 0.0; 5665e76c431SBarry Smith if (a == 0.0) { 5675e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5685e76c431SBarry Smith } else { 5695e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5705e76c431SBarry Smith } 5715e76c431SBarry Smith lambdaprev = lambda; 5725e76c431SBarry Smith gnormprev = *gnorm; 573329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 574329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5755e76c431SBarry Smith else lambda = lambdatemp; 576393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 577ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 578aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 579ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 580393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5815e42470aSBarry Smith #else 582ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5835e42470aSBarry Smith #endif 58478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 585cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5865d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 587393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 588b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5892715565aSLois Curfman McInnes goto theend1; 5902b022350SLois Curfman McInnes } else { 591b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5925e76c431SBarry Smith } 5935e76c431SBarry Smith count++; 5945e76c431SBarry Smith } 595ad922ac9SBarry Smith theend1: 5968f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5978f99978dSLois Curfman McInnes if (neP->CheckStep) { 5988f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5998f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 6008f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 6018f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6028f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 6038f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6048f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6058f99978dSLois Curfman McInnes } 6068f99978dSLois Curfman McInnes } 607b0a32e0cSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6083a40ed3dSBarry Smith PetscFunctionReturn(0); 6095e76c431SBarry Smith } 61004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6114a2ae208SSatish Balay #undef __FUNCT__ 6124a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch" 6134b828684SBarry Smith /*@C 614f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6155e76c431SBarry Smith 616c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 617c7afd0dbSLois Curfman McInnes 6185e42470aSBarry Smith Input Parameters: 619c7afd0dbSLois Curfman McInnes + snes - the SNES context 620af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6215e42470aSBarry Smith . x - current iterate 6225e42470aSBarry Smith . f - residual evaluated at x 6235e42470aSBarry Smith . y - search direction (contains new iterate on output) 6245e42470aSBarry Smith . w - work vector 625c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6265e42470aSBarry Smith 627c4a48953SLois Curfman McInnes Output Parameters: 628c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6295e42470aSBarry Smith . y - new iterate (contains search direction on input) 6305e42470aSBarry Smith . gnorm - 2-norm of g 6315e42470aSBarry Smith . ynorm - 2-norm of search length 632c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 633fee21e36SBarry Smith 634c4a48953SLois Curfman McInnes Options Database Key: 635c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 636c4a48953SLois Curfman McInnes 6375e42470aSBarry Smith Notes: 6386831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 6395e42470aSBarry Smith 64036851e7fSLois Curfman McInnes Level: advanced 64136851e7fSLois Curfman McInnes 642f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 643f59ffdedSLois Curfman McInnes 644af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6455e42470aSBarry Smith @*/ 6465d1a10b1SBarry 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) 6475e76c431SBarry Smith { 648406556e6SLois Curfman McInnes /* 649406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 650406556e6SLois Curfman McInnes minimization problem: 651406556e6SLois Curfman McInnes min z(x): R^n -> R, 652406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 653406556e6SLois Curfman McInnes */ 654329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 655aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 656ea709b57SSatish Balay PetscScalar cinitslope,clambda; 6576b5873e3SBarry Smith #endif 6585e42470aSBarry Smith int ierr,count; 6596831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 660ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 6618f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6625e76c431SBarry Smith 6633a40ed3dSBarry Smith PetscFunctionBegin; 664b0a32e0cSBarry Smith ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 665761aaf1bSLois Curfman McInnes *flag = 0; 6665e42470aSBarry Smith alpha = neP->alpha; 6675e42470aSBarry Smith maxstep = neP->maxstep; 6685e42470aSBarry Smith steptol = neP->steptol; 6695e76c431SBarry Smith 6703505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 671b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 672b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 673b37302c6SLois Curfman McInnes *gnorm = fnorm; 674b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 675b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 676ad922ac9SBarry Smith goto theend2; 67794a9d846SBarry Smith } 6785e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6795e42470aSBarry Smith scale = maxstep/(*ynorm); 680393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6815e42470aSBarry Smith *ynorm = maxstep; 6825e76c431SBarry Smith } 6835e42470aSBarry Smith minlambda = steptol/(*ynorm); 684a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 685aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 686a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 687329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6885e42470aSBarry Smith #else 689a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6905e42470aSBarry Smith #endif 6915e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6925e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6935e42470aSBarry Smith 694393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6955334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 69678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 697cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 6985d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 699393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 700b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 701ad922ac9SBarry Smith goto theend2; 7025e42470aSBarry Smith } 7035e42470aSBarry Smith 7045e42470aSBarry Smith /* Fit points with quadratic */ 705313b5538SBarry Smith lambda = 1.0; 7065e42470aSBarry Smith count = 1; 7075e42470aSBarry Smith while (1) { 7085e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 709b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 710b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 711393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 712761aaf1bSLois Curfman McInnes *flag = -1; break; 7135e42470aSBarry Smith } 714a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 715329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 716329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 717329f5518SBarry Smith else lambda = lambdatemp; 718393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7193505fcc1SBarry Smith lambdaneg = -lambda; 720aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7213505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7225e42470aSBarry Smith #else 7233505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7245e42470aSBarry Smith #endif 72578b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 726cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7275d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 728393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 729b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 73006259719SBarry Smith break; 7315e42470aSBarry Smith } 7325e42470aSBarry Smith count++; 7335e42470aSBarry Smith } 734ad922ac9SBarry Smith theend2: 7358f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7368f99978dSLois Curfman McInnes if (neP->CheckStep) { 7378f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7388f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7398f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7408f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7418f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7428f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7438f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7448f99978dSLois Curfman McInnes } 7458f99978dSLois Curfman McInnes } 746b0a32e0cSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7473a40ed3dSBarry Smith PetscFunctionReturn(0); 7485e76c431SBarry Smith } 74904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7504a2ae208SSatish Balay #undef __FUNCT__ 7514a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 752c9e6a524SLois Curfman McInnes /*@C 75377c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7546831982aSBarry Smith by the method SNESEQLS. 7555e76c431SBarry Smith 7565e76c431SBarry Smith Input Parameters: 757c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 758af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 759c7afd0dbSLois Curfman McInnes - func - pointer to int function 7605e76c431SBarry Smith 761fee21e36SBarry Smith Collective on SNES 762fee21e36SBarry Smith 763c4a48953SLois Curfman McInnes Available Routines: 764c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 765f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 766f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 767af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7685e76c431SBarry Smith 769c4a48953SLois Curfman McInnes Options Database Keys: 770af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 771c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 772c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 773c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 774c4a48953SLois Curfman McInnes 7755e76c431SBarry Smith Calling sequence of func: 776bd208895SLois Curfman McInnes .vb 777af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 778329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 779bd208895SLois Curfman McInnes .ve 7805e76c431SBarry Smith 7815e76c431SBarry Smith Input parameters for func: 782c7afd0dbSLois Curfman McInnes + snes - nonlinear context 783af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7845e76c431SBarry Smith . x - current iterate 7855e76c431SBarry Smith . f - residual evaluated at x 7865e76c431SBarry Smith . y - search direction (contains new iterate on output) 7875e76c431SBarry Smith . w - work vector 788c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7895e76c431SBarry Smith 7905e76c431SBarry Smith Output parameters for func: 791c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7925e76c431SBarry Smith . y - new iterate (contains search direction on input) 7935e76c431SBarry Smith . gnorm - 2-norm of g 7945e76c431SBarry Smith . ynorm - 2-norm of search length 795c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 796761aaf1bSLois Curfman McInnes on failure. 797f59ffdedSLois Curfman McInnes 79836851e7fSLois Curfman McInnes Level: advanced 79936851e7fSLois Curfman McInnes 800f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 801f59ffdedSLois Curfman McInnes 8025d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 8035d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 8045e76c431SBarry Smith @*/ 805329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 8065e76c431SBarry Smith { 807329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 80882bf6240SBarry Smith 8093a40ed3dSBarry Smith PetscFunctionBegin; 810b9617806SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)())&f);CHKERRQ(ierr); 81182bf6240SBarry Smith if (f) { 812af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 81382bf6240SBarry Smith } 8143a40ed3dSBarry Smith PetscFunctionReturn(0); 8155e76c431SBarry Smith } 81604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 817fb2e594dSBarry Smith EXTERN_C_BEGIN 8184a2ae208SSatish Balay #undef __FUNCT__ 8194a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS" 820af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 82187828ca2SBarry Smith PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 82282bf6240SBarry Smith { 82382bf6240SBarry Smith PetscFunctionBegin; 8246831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 8256831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 82682bf6240SBarry Smith PetscFunctionReturn(0); 82782bf6240SBarry Smith } 828fb2e594dSBarry Smith EXTERN_C_END 82904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8304a2ae208SSatish Balay #undef __FUNCT__ 8314a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck" 832c8dd0c0dSLois Curfman McInnes /*@C 833530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 8346831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 835c8dd0c0dSLois Curfman McInnes 836c8dd0c0dSLois Curfman McInnes Input Parameters: 837c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 838c8dd0c0dSLois Curfman McInnes . func - pointer to int function 839c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 840c8dd0c0dSLois Curfman McInnes 841c8dd0c0dSLois Curfman McInnes Collective on SNES 842c8dd0c0dSLois Curfman McInnes 843c8dd0c0dSLois Curfman McInnes Calling sequence of func: 844c8dd0c0dSLois Curfman McInnes .vb 845b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 846c8dd0c0dSLois Curfman McInnes .ve 847b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 848b1ae27eaSLois Curfman McInnes on failure. 849c8dd0c0dSLois Curfman McInnes 850c8dd0c0dSLois Curfman McInnes Input parameters for func: 851c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 852c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8538f99978dSLois Curfman McInnes - x - current candidate iterate 854c8dd0c0dSLois Curfman McInnes 855c8dd0c0dSLois Curfman McInnes Output parameters for func: 856c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 857c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 858c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 859c8dd0c0dSLois Curfman McInnes 860c8dd0c0dSLois Curfman McInnes Level: advanced 861c8dd0c0dSLois Curfman McInnes 8628f99978dSLois Curfman McInnes Notes: 863b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 864b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 865b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 866b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 867ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 8688f99978dSLois Curfman McInnes 869b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 870b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 871b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 872ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 873ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 874ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 875ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 876b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8778f99978dSLois Curfman McInnes 878c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 879c8dd0c0dSLois Curfman McInnes 880c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 881c8dd0c0dSLois Curfman McInnes @*/ 8828f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 883c8dd0c0dSLois Curfman McInnes { 8848f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 885c8dd0c0dSLois Curfman McInnes 886c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 887b9617806SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)())&f);CHKERRQ(ierr); 888c8dd0c0dSLois Curfman McInnes if (f) { 889c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 890c8dd0c0dSLois Curfman McInnes } 891c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 892c8dd0c0dSLois Curfman McInnes } 893c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 894c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 8954a2ae208SSatish Balay #undef __FUNCT__ 8964a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS" 8978f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 898c8dd0c0dSLois Curfman McInnes { 899c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 9006831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 9016831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 902c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 903c8dd0c0dSLois Curfman McInnes } 904c8dd0c0dSLois Curfman McInnes EXTERN_C_END 905c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 90604d965bbSLois Curfman McInnes /* 907*329e7e40SMatthew Knepley SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method. 908*329e7e40SMatthew Knepley 909*329e7e40SMatthew Knepley Input Parameter: 910*329e7e40SMatthew Knepley . snes - the SNES context 911*329e7e40SMatthew Knepley 912*329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 913*329e7e40SMatthew Knepley */ 914*329e7e40SMatthew Knepley #undef __FUNCT__ 915*329e7e40SMatthew Knepley #define __FUNCT__ "SNESPrintHelp_EQ_LS" 916*329e7e40SMatthew Knepley static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 917*329e7e40SMatthew Knepley { 918*329e7e40SMatthew Knepley SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 919*329e7e40SMatthew Knepley 920*329e7e40SMatthew Knepley PetscFunctionBegin; 921*329e7e40SMatthew Knepley (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n"); 922*329e7e40SMatthew Knepley (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 923*329e7e40SMatthew Knepley (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 924*329e7e40SMatthew Knepley (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 925*329e7e40SMatthew Knepley (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol); 926*329e7e40SMatthew Knepley PetscFunctionReturn(0); 927*329e7e40SMatthew Knepley } 928*329e7e40SMatthew Knepley 929*329e7e40SMatthew Knepley /* 9306831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 93104d965bbSLois Curfman McInnes 93204d965bbSLois Curfman McInnes Input Parameters: 93304d965bbSLois Curfman McInnes . SNES - the SNES context 93404d965bbSLois Curfman McInnes . viewer - visualization context 93504d965bbSLois Curfman McInnes 93604d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 93704d965bbSLois Curfman McInnes */ 9384a2ae208SSatish Balay #undef __FUNCT__ 9394a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_LS" 940b0a32e0cSBarry Smith static int SNESView_EQ_LS(SNES snes,PetscViewer viewer) 941a935fc98SLois Curfman McInnes { 9426831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 94319bcc07fSBarry Smith char *cstr; 94451695f54SLois Curfman McInnes int ierr; 9456831982aSBarry Smith PetscTruth isascii; 946a935fc98SLois Curfman McInnes 9473a40ed3dSBarry Smith PetscFunctionBegin; 948b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 9490f5bd95cSBarry Smith if (isascii) { 95019bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 95119bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 95219bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 95319bcc07fSBarry Smith else cstr = "unknown"; 954b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 955b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9565cd90555SBarry Smith } else { 95729bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 95819bcc07fSBarry Smith } 9593a40ed3dSBarry Smith PetscFunctionReturn(0); 960a935fc98SLois Curfman McInnes } 96104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 96204d965bbSLois Curfman McInnes /* 9636831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 96404d965bbSLois Curfman McInnes 96504d965bbSLois Curfman McInnes Input Parameter: 96604d965bbSLois Curfman McInnes . snes - the SNES context 96704d965bbSLois Curfman McInnes 96804d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 96904d965bbSLois Curfman McInnes */ 9704a2ae208SSatish Balay #undef __FUNCT__ 9714a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_LS" 972f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 9735e42470aSBarry Smith { 9746831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 975186905e3SBarry Smith char ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"}; 976f1af5d2fSBarry Smith int ierr; 977f1af5d2fSBarry Smith PetscTruth flg; 9785e42470aSBarry Smith 9793a40ed3dSBarry Smith PetscFunctionBegin; 980b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 98187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 98287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 98387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 984186905e3SBarry Smith 985b0a32e0cSBarry Smith ierr = PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);CHKERRQ(ierr); 98625cdf11fSBarry Smith if (flg) { 987c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9880f5bd95cSBarry Smith 989186905e3SBarry Smith ierr = PetscStrcmp(ver,lses[0],&isbasic);CHKERRQ(ierr); 990186905e3SBarry Smith ierr = PetscStrcmp(ver,lses[1],&isnonorms);CHKERRQ(ierr); 991186905e3SBarry Smith ierr = PetscStrcmp(ver,lses[2],&isquad);CHKERRQ(ierr); 992186905e3SBarry Smith ierr = PetscStrcmp(ver,lses[3],&iscubic);CHKERRQ(ierr); 9930f5bd95cSBarry Smith 9940f5bd95cSBarry Smith if (isbasic) { 995af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9960f5bd95cSBarry Smith } else if (isnonorms) { 997af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9980f5bd95cSBarry Smith } else if (isquad) { 999af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 10000f5bd95cSBarry Smith } else if (iscubic) { 1001af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 10025e42470aSBarry Smith } 100329bbc08cSBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");} 10045e42470aSBarry Smith } 1005b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10063a40ed3dSBarry Smith PetscFunctionReturn(0); 10075e42470aSBarry Smith } 100804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 100904d965bbSLois Curfman McInnes /* 10106831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 10116831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 101204d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 101304d965bbSLois Curfman McInnes 101404d965bbSLois Curfman McInnes 101504d965bbSLois Curfman McInnes Input Parameter: 101604d965bbSLois Curfman McInnes . snes - the SNES context 101704d965bbSLois Curfman McInnes 101804d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 101904d965bbSLois Curfman McInnes */ 1020fb2e594dSBarry Smith EXTERN_C_BEGIN 10214a2ae208SSatish Balay #undef __FUNCT__ 10224a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_LS" 1023f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 10245e42470aSBarry Smith { 102582bf6240SBarry Smith int ierr; 10266831982aSBarry Smith SNES_EQ_LS *neP; 10275e42470aSBarry Smith 10283a40ed3dSBarry Smith PetscFunctionBegin; 1029a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 103029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1031a8c6a408SBarry Smith } 103282bf6240SBarry Smith 1033f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 1034f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 1035f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 1036f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 1037*329e7e40SMatthew Knepley snes->printhelp = SNESPrintHelp_EQ_LS; 1038f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 1039f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 10405baf8537SBarry Smith snes->nwork = 0; 10415e42470aSBarry Smith 1042b0a32e0cSBarry Smith ierr = PetscNew(SNES_EQ_LS,&neP);CHKERRQ(ierr); 1043b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 10445e42470aSBarry Smith snes->data = (void*)neP; 10455e42470aSBarry Smith neP->alpha = 1.e-4; 10465e42470aSBarry Smith neP->maxstep = 1.e8; 10475e42470aSBarry Smith neP->steptol = 1.e-12; 10485e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1049c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1050c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1051c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 105282bf6240SBarry Smith 10535d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10545d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 105582bf6240SBarry Smith 10563a40ed3dSBarry Smith PetscFunctionReturn(0); 10575e42470aSBarry Smith } 1058fb2e594dSBarry Smith EXTERN_C_END 105982bf6240SBarry Smith 106082bf6240SBarry Smith 106182bf6240SBarry Smith 1062