1*74637425SBarry Smith /*$Id: ls.c,v 1.159 2000/07/28 14:46:38 bsmith Exp bsmith $*/ 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 */ 118a5d9ee4SBarry Smith #undef __FUNC__ 128a5d9ee4SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"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); 265d1a10b1SBarry Smith PLogInfo(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; 3036669109SBarry Smith Scalar 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); 405d1a10b1SBarry Smith PLogInfo(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 46*74637425SBarry Smith /* 47*74637425SBarry Smith Checks if J^T(F - AX) = 0 48*74637425SBarry Smith */ 49*74637425SBarry Smith #undef __FUNC__ 50*74637425SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckResidual_Private"></a>*/"SNESLSCheckResidual_Private" 51*74637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 52*74637425SBarry Smith { 53*74637425SBarry Smith PetscReal a1,a2; 54*74637425SBarry Smith int ierr; 55*74637425SBarry Smith PetscTruth hastranspose; 56*74637425SBarry Smith Scalar mone = -1.0; 57*74637425SBarry Smith 58*74637425SBarry Smith PetscFunctionBegin; 59*74637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 60*74637425SBarry Smith if (hastranspose) { 61*74637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 62*74637425SBarry Smith ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr); 63*74637425SBarry Smith 64*74637425SBarry Smith /* Compute || J^T W|| */ 65*74637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 66*74637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 67*74637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 68*74637425SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: || J^T(F - Ax)|| %g near zero implies inconsistent rhs\n",a2/a1); 69*74637425SBarry Smith } 70*74637425SBarry Smith PetscFunctionReturn(0); 71*74637425SBarry Smith } 72*74637425SBarry Smith 7304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7404d965bbSLois Curfman McInnes 7504d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7604d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 7704d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7804d965bbSLois Curfman McInnes respectively. 7904d965bbSLois Curfman McInnes 80fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8104d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8204d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 83fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8404d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8504d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 8604d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 8704d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 8804d965bbSLois Curfman McInnes These routines are actually called via the common user interface 8904d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9004d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9104d965bbSLois Curfman McInnes for all nonlinear solvers. 9204d965bbSLois Curfman McInnes 9304d965bbSLois Curfman McInnes Another key routine is: 9404d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9504d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9604d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9704d965bbSLois Curfman McInnes SNESSolve() if necessary. 9804d965bbSLois Curfman McInnes 9904d965bbSLois Curfman McInnes Additional basic routines are: 10004d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 10104d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10204d965bbSLois Curfman McInnes have actually been used. 10304d965bbSLois Curfman McInnes These are called by application codes via the interface routines 10404d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 10504d965bbSLois Curfman McInnes 10604d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10704d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10804d965bbSLois Curfman McInnes above description applies to these categories also. 10904d965bbSLois Curfman McInnes 11004d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1115e42470aSBarry Smith /* 11204d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 11304d965bbSLois Curfman McInnes method with a line search. 1145e76c431SBarry Smith 11504d965bbSLois Curfman McInnes Input Parameters: 11604d965bbSLois Curfman McInnes . snes - the SNES context 1175e76c431SBarry Smith 11804d965bbSLois Curfman McInnes Output Parameter: 119c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12004d965bbSLois Curfman McInnes 12104d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1225e76c431SBarry Smith 1235e76c431SBarry Smith Notes: 1245e76c431SBarry Smith This implements essentially a truncated Newton method with a 1255e76c431SBarry Smith line search. By default a cubic backtracking line search 1265e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1275e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 128393d2d9aSLois Curfman McInnes and Schnabel. 1295e76c431SBarry Smith */ 1305615d1e5SSatish Balay #undef __FUNC__ 131b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 132f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 1335e76c431SBarry Smith { 1346831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 13584c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 136112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 137329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1385e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 1395e76c431SBarry Smith 1403a40ed3dSBarry Smith PetscFunctionBegin; 141184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 142184914b5SBarry Smith 1435e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1445e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1465e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1475e42470aSBarry Smith G = snes->work[1]; 1485e42470aSBarry Smith W = snes->work[2]; 1495e76c431SBarry Smith 1504c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1514c49b128SBarry Smith snes->iter = 0; 1524c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1535334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 154cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1550f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1565e42470aSBarry Smith snes->norm = fnorm; 1570f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15884c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 15994a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1605e76c431SBarry Smith 161184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16294a9d846SBarry Smith 163d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 164d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 165d034289bSBarry Smith 1665e76c431SBarry Smith for (i=0; i<maxits; i++) { 1675e76c431SBarry Smith 168ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1695334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1705334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 17178b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 172*74637425SBarry Smith 173*74637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 174*74637425SBarry Smith 17590cbd584SBarry Smith /* should check what happened to the linear solve? */ 1763505fcc1SBarry Smith snes->linear_its += lits; 177af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 178ea4d3ed3SLois Curfman McInnes 179ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 180ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 181ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 182ea4d3ed3SLois Curfman McInnes */ 18381b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 184af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 185af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 1863505fcc1SBarry Smith if (lsfail) { 1878a5d9ee4SBarry Smith PetscTruth ismin; 1883505fcc1SBarry Smith snes->nfailures++; 1893505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 1908a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 1918a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1923505fcc1SBarry Smith break; 1933505fcc1SBarry Smith } 1945e76c431SBarry Smith 19539e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 19639e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1975e76c431SBarry Smith fnorm = gnorm; 1985e76c431SBarry Smith 1990f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 20084c569c9SBarry Smith snes->iter = i+1; 2015e42470aSBarry Smith snes->norm = fnorm; 2020f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 20384c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 20494a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 2055e76c431SBarry Smith 2065e76c431SBarry Smith /* Test for convergence */ 20729e0b56fSBarry Smith if (snes->converged) { 20829e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2093505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2103505fcc1SBarry Smith if (snes->reason) { 21116c95cacSBarry Smith break; 21216c95cacSBarry Smith } 21316c95cacSBarry Smith } 21429e0b56fSBarry Smith } 21539e2f89bSBarry Smith if (X != snes->vec_sol) { 216393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 217e15f7bb5SBarry Smith } 218e15f7bb5SBarry Smith if (F != snes->vec_func) { 219e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 220e15f7bb5SBarry Smith } 22139e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 22239e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 22352392280SLois Curfman McInnes if (i == maxits) { 224981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 22552392280SLois Curfman McInnes i--; 2263505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 22752392280SLois Curfman McInnes } 2280f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2290f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2305e42470aSBarry Smith *outits = i+1; 2313a40ed3dSBarry Smith PetscFunctionReturn(0); 2325e76c431SBarry Smith } 23304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 23404d965bbSLois Curfman McInnes /* 23504d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 2366831982aSBarry Smith of the SNESEQLS nonlinear solver. 23704d965bbSLois Curfman McInnes 23804d965bbSLois Curfman McInnes Input Parameter: 23904d965bbSLois Curfman McInnes . snes - the SNES context 24004d965bbSLois Curfman McInnes . x - the solution vector 24104d965bbSLois Curfman McInnes 24204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 24304d965bbSLois Curfman McInnes 24404d965bbSLois Curfman McInnes Notes: 24504d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 24604d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 24704d965bbSLois Curfman McInnes the call to SNESSolve(). 24804d965bbSLois Curfman McInnes */ 2495615d1e5SSatish Balay #undef __FUNC__ 250b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 251f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 2525e76c431SBarry Smith { 2535e42470aSBarry Smith int ierr; 2543a40ed3dSBarry Smith 2553a40ed3dSBarry Smith PetscFunctionBegin; 25681b6cf68SLois Curfman McInnes snes->nwork = 4; 257d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 2585e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 25981b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2603a40ed3dSBarry Smith PetscFunctionReturn(0); 2615e76c431SBarry Smith } 26204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 26304d965bbSLois Curfman McInnes /* 2646831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 26504d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 26604d965bbSLois Curfman McInnes 26704d965bbSLois Curfman McInnes Input Parameter: 26804d965bbSLois Curfman McInnes . snes - the SNES context 26904d965bbSLois Curfman McInnes 270de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 27104d965bbSLois Curfman McInnes */ 2725615d1e5SSatish Balay #undef __FUNC__ 273b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 274e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2755e76c431SBarry Smith { 276393d2d9aSLois Curfman McInnes int ierr; 2773a40ed3dSBarry Smith 2783a40ed3dSBarry Smith PetscFunctionBegin; 2795baf8537SBarry Smith if (snes->nwork) { 2804b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 28121c89e3eSBarry Smith } 2825c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2845e76c431SBarry Smith } 28504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2865615d1e5SSatish Balay #undef __FUNC__ 287b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 28882bf6240SBarry Smith 2894b828684SBarry Smith /*@C 2905e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2915e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2925e76c431SBarry Smith to serve as a template and is not recommended for general use. 2935e76c431SBarry Smith 294c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 295c7afd0dbSLois Curfman McInnes 2965e76c431SBarry Smith Input Parameters: 297c7afd0dbSLois Curfman McInnes + snes - nonlinear context 298af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2995e76c431SBarry Smith . x - current iterate 3005e76c431SBarry Smith . f - residual evaluated at x 3015e76c431SBarry Smith . y - search direction (contains new iterate on output) 3025e76c431SBarry Smith . w - work vector 303c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3045e76c431SBarry Smith 305c4a48953SLois Curfman McInnes Output Parameters: 306c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3075e76c431SBarry Smith . y - new iterate (contains search direction on input) 3085e76c431SBarry Smith . gnorm - 2-norm of g 3095e76c431SBarry Smith . ynorm - 2-norm of search length 310c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 311fee21e36SBarry Smith 312c4a48953SLois Curfman McInnes Options Database Key: 313c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 314c4a48953SLois Curfman McInnes 31536851e7fSLois Curfman McInnes Level: advanced 31636851e7fSLois Curfman McInnes 31728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 31828ae5a21SLois Curfman McInnes 319f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 320af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3215e76c431SBarry Smith @*/ 3225d1a10b1SBarry 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) 3235e76c431SBarry Smith { 3245e42470aSBarry Smith int ierr; 3255334005bSBarry Smith Scalar mone = -1.0; 3266831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3278f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3285334005bSBarry Smith 3293a40ed3dSBarry Smith PetscFunctionBegin; 330761aaf1bSLois Curfman McInnes *flag = 0; 3317857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 332a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 333ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3348f99978dSLois Curfman McInnes if (neP->CheckStep) { 3358f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3368f99978dSLois Curfman McInnes } 337ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 338a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 3397857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3403a40ed3dSBarry Smith PetscFunctionReturn(0); 3415e76c431SBarry Smith } 34204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3435615d1e5SSatish Balay #undef __FUNC__ 344b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 34582bf6240SBarry Smith 34629e0b56fSBarry Smith /*@C 34729e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 34829e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 34929e0b56fSBarry Smith even compute the norm of the function or search direction; this 35029e0b56fSBarry Smith is intended only when you know the full step is fine and are 35129e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 352c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 353c7afd0dbSLois Curfman McInnes 354c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 35529e0b56fSBarry Smith 35629e0b56fSBarry Smith Input Parameters: 357c7afd0dbSLois Curfman McInnes + snes - nonlinear context 358af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 35929e0b56fSBarry Smith . x - current iterate 36029e0b56fSBarry Smith . f - residual evaluated at x 36129e0b56fSBarry Smith . y - search direction (contains new iterate on output) 36229e0b56fSBarry Smith . w - work vector 363c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 36429e0b56fSBarry Smith 36529e0b56fSBarry Smith Output Parameters: 366c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 36729e0b56fSBarry Smith . gnorm - not changed 36829e0b56fSBarry Smith . ynorm - not changed 369c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 370fee21e36SBarry Smith 37129e0b56fSBarry Smith Options Database Key: 372c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 37329e0b56fSBarry Smith 3748cbba510SBarry Smith Notes: 375ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 376ea56f5baSLois Curfman McInnes either the options 377ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 378ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 379329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 380329f5518SBarry Smith otherwise, the SNES solver will generate an error. 381329f5518SBarry Smith 382329f5518SBarry Smith During the final iteration this will not evaluate the function at 383329f5518SBarry Smith the solution point. This is to save a function evaluation while 384329f5518SBarry Smith using pseudo-timestepping. 3858cbba510SBarry Smith 386ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 387ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 388ea56f5baSLois Curfman McInnes correct, since they are not computed. 389ea56f5baSLois Curfman McInnes 390ea56f5baSLois Curfman McInnes Level: advanced 3918cbba510SBarry Smith 39229e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 39329e0b56fSBarry Smith 39429e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 39529e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 39629e0b56fSBarry Smith @*/ 3975d1a10b1SBarry 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) 39829e0b56fSBarry Smith { 39929e0b56fSBarry Smith int ierr; 40029e0b56fSBarry Smith Scalar mone = -1.0; 4016831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4028f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 40329e0b56fSBarry Smith 4043a40ed3dSBarry Smith PetscFunctionBegin; 4058cbba510SBarry Smith *flag = 0; 40629e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 40729e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4088f99978dSLois Curfman McInnes if (neP->CheckStep) { 4098f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4108f99978dSLois Curfman McInnes } 411329f5518SBarry Smith 412329f5518SBarry Smith /* don't evaluate function the last time through */ 413329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 41429e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 415329f5518SBarry Smith } 41629e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4173a40ed3dSBarry Smith PetscFunctionReturn(0); 41829e0b56fSBarry Smith } 41904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 42029e0b56fSBarry Smith #undef __FUNC__ 421b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 4224b828684SBarry Smith /*@C 423f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4245e76c431SBarry Smith 425c7afd0dbSLois Curfman McInnes Collective on SNES 426c7afd0dbSLois Curfman McInnes 4275e76c431SBarry Smith Input Parameters: 428c7afd0dbSLois Curfman McInnes + snes - nonlinear context 429af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4305e76c431SBarry Smith . x - current iterate 4315e76c431SBarry Smith . f - residual evaluated at x 4325e76c431SBarry Smith . y - search direction (contains new iterate on output) 4335e76c431SBarry Smith . w - work vector 434c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4355e76c431SBarry Smith 436393d2d9aSLois Curfman McInnes Output Parameters: 437c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4385e76c431SBarry Smith . y - new iterate (contains search direction on input) 4395e76c431SBarry Smith . gnorm - 2-norm of g 4405e76c431SBarry Smith . ynorm - 2-norm of search length 441c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 442fee21e36SBarry Smith 443c4a48953SLois Curfman McInnes Options Database Key: 444c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 445c4a48953SLois Curfman McInnes 4465e76c431SBarry Smith Notes: 4475e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4485e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4495e76c431SBarry Smith 45036851e7fSLois Curfman McInnes Level: advanced 45136851e7fSLois Curfman McInnes 45228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 45328ae5a21SLois Curfman McInnes 454af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4555e76c431SBarry Smith @*/ 4565d1a10b1SBarry 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) 4575e76c431SBarry Smith { 458406556e6SLois Curfman McInnes /* 459406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 460406556e6SLois Curfman McInnes minimization problem: 461406556e6SLois Curfman McInnes min z(x): R^n -> R, 462406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 463406556e6SLois Curfman McInnes */ 464406556e6SLois Curfman McInnes 465329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 466329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 467aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4685e42470aSBarry Smith Scalar cinitslope,clambda; 4696b5873e3SBarry Smith #endif 4705e42470aSBarry Smith int ierr,count; 4716831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4725334005bSBarry Smith Scalar mone = -1.0,scale; 4738f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4745e76c431SBarry Smith 4753a40ed3dSBarry Smith PetscFunctionBegin; 4767857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 477761aaf1bSLois Curfman McInnes *flag = 0; 4785e76c431SBarry Smith alpha = neP->alpha; 4795e76c431SBarry Smith maxstep = neP->maxstep; 4805e76c431SBarry Smith steptol = neP->steptol; 4815e76c431SBarry Smith 482cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 483a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 484a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 485a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 486a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 487a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 488ad922ac9SBarry Smith goto theend1; 48994a9d846SBarry Smith } 4905e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4915e42470aSBarry Smith scale = maxstep/(*ynorm); 492aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 49390cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 4946b5873e3SBarry Smith #else 49590cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 4966b5873e3SBarry Smith #endif 497393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4985e76c431SBarry Smith *ynorm = maxstep; 4995e76c431SBarry Smith } 5005e76c431SBarry Smith minlambda = steptol/(*ynorm); 501a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 502aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 503a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 504329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5055e42470aSBarry Smith #else 506a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5075e42470aSBarry Smith #endif 5085e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5095e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5105e76c431SBarry Smith 511393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5125334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 51378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 514313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5155d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 516393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 51794a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 518ad922ac9SBarry Smith goto theend1; 5195e76c431SBarry Smith } 5205e76c431SBarry Smith 5215e76c431SBarry Smith /* Fit points with quadratic */ 522313b5538SBarry Smith lambda = 1.0; 523a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5245e76c431SBarry Smith lambdaprev = lambda; 5255e76c431SBarry Smith gnormprev = *gnorm; 526329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 527ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 528ddd12767SBarry Smith else lambda = lambdatemp; 529393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 530ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 531aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 532ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5335e42470aSBarry Smith #else 534ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5355e42470aSBarry Smith #endif 53678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 537cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5385d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 539393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 540ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 541ad922ac9SBarry Smith goto theend1; 5425e76c431SBarry Smith } 5435e76c431SBarry Smith 5445e76c431SBarry Smith /* Fit points with cubic */ 5455e76c431SBarry Smith count = 1; 5465e76c431SBarry Smith while (1) { 5475e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 5482b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 54990cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 550393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 551761aaf1bSLois Curfman McInnes *flag = -1; break; 5525e76c431SBarry Smith } 5535d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 5545d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 555ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5562b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5575e76c431SBarry Smith d = b*b - 3*a*initslope; 5585e76c431SBarry Smith if (d < 0.0) d = 0.0; 5595e76c431SBarry Smith if (a == 0.0) { 5605e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5615e76c431SBarry Smith } else { 5625e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5635e76c431SBarry Smith } 5645e76c431SBarry Smith lambdaprev = lambda; 5655e76c431SBarry Smith gnormprev = *gnorm; 566329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 567329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5685e76c431SBarry Smith else lambda = lambdatemp; 569393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 570ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 571aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 572ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 573393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5745e42470aSBarry Smith #else 575ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5765e42470aSBarry Smith #endif 57778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 578cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 5795d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 580393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 581ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5822715565aSLois Curfman McInnes goto theend1; 5832b022350SLois Curfman McInnes } else { 5842b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5855e76c431SBarry Smith } 5865e76c431SBarry Smith count++; 5875e76c431SBarry Smith } 588ad922ac9SBarry Smith theend1: 5898f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5908f99978dSLois Curfman McInnes if (neP->CheckStep) { 5918f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5928f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5938f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5948f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5958f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5968f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5978f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5988f99978dSLois Curfman McInnes } 5998f99978dSLois Curfman McInnes } 6007857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 6013a40ed3dSBarry Smith PetscFunctionReturn(0); 6025e76c431SBarry Smith } 60304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6045615d1e5SSatish Balay #undef __FUNC__ 605b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 6064b828684SBarry Smith /*@C 607f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6085e76c431SBarry Smith 609c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 610c7afd0dbSLois Curfman McInnes 6115e42470aSBarry Smith Input Parameters: 612c7afd0dbSLois Curfman McInnes + snes - the SNES context 613af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6145e42470aSBarry Smith . x - current iterate 6155e42470aSBarry Smith . f - residual evaluated at x 6165e42470aSBarry Smith . y - search direction (contains new iterate on output) 6175e42470aSBarry Smith . w - work vector 618c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6195e42470aSBarry Smith 620c4a48953SLois Curfman McInnes Output Parameters: 621c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6225e42470aSBarry Smith . y - new iterate (contains search direction on input) 6235e42470aSBarry Smith . gnorm - 2-norm of g 6245e42470aSBarry Smith . ynorm - 2-norm of search length 625c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 626fee21e36SBarry Smith 627c4a48953SLois Curfman McInnes Options Database Key: 628c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 629c4a48953SLois Curfman McInnes 6305e42470aSBarry Smith Notes: 6316831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 6325e42470aSBarry Smith 63336851e7fSLois Curfman McInnes Level: advanced 63436851e7fSLois Curfman McInnes 635f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 636f59ffdedSLois Curfman McInnes 637af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6385e42470aSBarry Smith @*/ 6395d1a10b1SBarry 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) 6405e76c431SBarry Smith { 641406556e6SLois Curfman McInnes /* 642406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 643406556e6SLois Curfman McInnes minimization problem: 644406556e6SLois Curfman McInnes min z(x): R^n -> R, 645406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 646406556e6SLois Curfman McInnes */ 647329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 648aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6495e42470aSBarry Smith Scalar cinitslope,clambda; 6506b5873e3SBarry Smith #endif 6515e42470aSBarry Smith int ierr,count; 6526831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 6535334005bSBarry Smith Scalar mone = -1.0,scale; 6548f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6555e76c431SBarry Smith 6563a40ed3dSBarry Smith PetscFunctionBegin; 6577857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 658761aaf1bSLois Curfman McInnes *flag = 0; 6595e42470aSBarry Smith alpha = neP->alpha; 6605e42470aSBarry Smith maxstep = neP->maxstep; 6615e42470aSBarry Smith steptol = neP->steptol; 6625e76c431SBarry Smith 6633505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 664b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 66594a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 666b37302c6SLois Curfman McInnes *gnorm = fnorm; 667b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 668b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 669ad922ac9SBarry Smith goto theend2; 67094a9d846SBarry Smith } 6715e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6725e42470aSBarry Smith scale = maxstep/(*ynorm); 673393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6745e42470aSBarry Smith *ynorm = maxstep; 6755e76c431SBarry Smith } 6765e42470aSBarry Smith minlambda = steptol/(*ynorm); 677a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 678aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 679a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 680329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6815e42470aSBarry Smith #else 682a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6835e42470aSBarry Smith #endif 6845e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6855e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6865e42470aSBarry Smith 687393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6885334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 68978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 690cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 6915d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 692393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 69394a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 694ad922ac9SBarry Smith goto theend2; 6955e42470aSBarry Smith } 6965e42470aSBarry Smith 6975e42470aSBarry Smith /* Fit points with quadratic */ 698313b5538SBarry Smith lambda = 1.0; 6995e42470aSBarry Smith count = 1; 7005e42470aSBarry Smith while (1) { 7015e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 702981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 7035d1a10b1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 704393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 705761aaf1bSLois Curfman McInnes *flag = -1; break; 7065e42470aSBarry Smith } 707a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 708329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 709329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 710329f5518SBarry Smith else lambda = lambdatemp; 711393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 7123505fcc1SBarry Smith lambdaneg = -lambda; 713aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 7143505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 7155e42470aSBarry Smith #else 7163505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 7175e42470aSBarry Smith #endif 71878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 719cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 7205d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 721393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 722981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 72306259719SBarry Smith break; 7245e42470aSBarry Smith } 7255e42470aSBarry Smith count++; 7265e42470aSBarry Smith } 727ad922ac9SBarry Smith theend2: 7288f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 7298f99978dSLois Curfman McInnes if (neP->CheckStep) { 7308f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7318f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7328f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7338f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7348f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7358f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7368f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7378f99978dSLois Curfman McInnes } 7388f99978dSLois Curfman McInnes } 7397857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 7403a40ed3dSBarry Smith PetscFunctionReturn(0); 7415e76c431SBarry Smith } 74204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7435615d1e5SSatish Balay #undef __FUNC__ 744b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 745c9e6a524SLois Curfman McInnes /*@C 74677c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7476831982aSBarry Smith by the method SNESEQLS. 7485e76c431SBarry Smith 7495e76c431SBarry Smith Input Parameters: 750c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 751af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 752c7afd0dbSLois Curfman McInnes - func - pointer to int function 7535e76c431SBarry Smith 754fee21e36SBarry Smith Collective on SNES 755fee21e36SBarry Smith 756c4a48953SLois Curfman McInnes Available Routines: 757c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 758f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 759f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 760af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7615e76c431SBarry Smith 762c4a48953SLois Curfman McInnes Options Database Keys: 763af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 764c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 765c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 766c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 767c4a48953SLois Curfman McInnes 7685e76c431SBarry Smith Calling sequence of func: 769bd208895SLois Curfman McInnes .vb 770af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 771329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 772bd208895SLois Curfman McInnes .ve 7735e76c431SBarry Smith 7745e76c431SBarry Smith Input parameters for func: 775c7afd0dbSLois Curfman McInnes + snes - nonlinear context 776af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7775e76c431SBarry Smith . x - current iterate 7785e76c431SBarry Smith . f - residual evaluated at x 7795e76c431SBarry Smith . y - search direction (contains new iterate on output) 7805e76c431SBarry Smith . w - work vector 781c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7825e76c431SBarry Smith 7835e76c431SBarry Smith Output parameters for func: 784c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7855e76c431SBarry Smith . y - new iterate (contains search direction on input) 7865e76c431SBarry Smith . gnorm - 2-norm of g 7875e76c431SBarry Smith . ynorm - 2-norm of search length 788c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 789761aaf1bSLois Curfman McInnes on failure. 790f59ffdedSLois Curfman McInnes 79136851e7fSLois Curfman McInnes Level: advanced 79236851e7fSLois Curfman McInnes 793f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 794f59ffdedSLois Curfman McInnes 7955d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 7965d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7975e76c431SBarry Smith @*/ 798329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7995e76c431SBarry Smith { 800329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 80182bf6240SBarry Smith 8023a40ed3dSBarry Smith PetscFunctionBegin; 80394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 80482bf6240SBarry Smith if (f) { 805af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 80682bf6240SBarry Smith } 8073a40ed3dSBarry Smith PetscFunctionReturn(0); 8085e76c431SBarry Smith } 80904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 810fb2e594dSBarry Smith EXTERN_C_BEGIN 81182bf6240SBarry Smith #undef __FUNC__ 812b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 813af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 814af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 81582bf6240SBarry Smith { 81682bf6240SBarry Smith PetscFunctionBegin; 8176831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 8186831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 81982bf6240SBarry Smith PetscFunctionReturn(0); 82082bf6240SBarry Smith } 821fb2e594dSBarry Smith EXTERN_C_END 82204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 823c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 824b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 825c8dd0c0dSLois Curfman McInnes /*@C 826530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 8276831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 828c8dd0c0dSLois Curfman McInnes 829c8dd0c0dSLois Curfman McInnes Input Parameters: 830c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 831c8dd0c0dSLois Curfman McInnes . func - pointer to int function 832c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 833c8dd0c0dSLois Curfman McInnes 834c8dd0c0dSLois Curfman McInnes Collective on SNES 835c8dd0c0dSLois Curfman McInnes 836c8dd0c0dSLois Curfman McInnes Calling sequence of func: 837c8dd0c0dSLois Curfman McInnes .vb 838b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 839c8dd0c0dSLois Curfman McInnes .ve 840b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 841b1ae27eaSLois Curfman McInnes on failure. 842c8dd0c0dSLois Curfman McInnes 843c8dd0c0dSLois Curfman McInnes Input parameters for func: 844c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 845c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8468f99978dSLois Curfman McInnes - x - current candidate iterate 847c8dd0c0dSLois Curfman McInnes 848c8dd0c0dSLois Curfman McInnes Output parameters for func: 849c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 850c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 851c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 852c8dd0c0dSLois Curfman McInnes 853c8dd0c0dSLois Curfman McInnes Level: advanced 854c8dd0c0dSLois Curfman McInnes 8558f99978dSLois Curfman McInnes Notes: 856b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 857b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 858b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 859b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 860ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 8618f99978dSLois Curfman McInnes 862b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 863b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 864b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 865ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 866ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 867ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 868ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 869b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8708f99978dSLois Curfman McInnes 871c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 872c8dd0c0dSLois Curfman McInnes 873c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 874c8dd0c0dSLois Curfman McInnes @*/ 8758f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 876c8dd0c0dSLois Curfman McInnes { 8778f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 878c8dd0c0dSLois Curfman McInnes 879c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 880c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 881c8dd0c0dSLois Curfman McInnes if (f) { 882c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 883c8dd0c0dSLois Curfman McInnes } 884c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 885c8dd0c0dSLois Curfman McInnes } 886c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 887c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 888c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 889b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8908f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 891c8dd0c0dSLois Curfman McInnes { 892c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8936831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8946831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 895c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 896c8dd0c0dSLois Curfman McInnes } 897c8dd0c0dSLois Curfman McInnes EXTERN_C_END 898c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 89904d965bbSLois Curfman McInnes /* 9006831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 90182bf6240SBarry Smith 90204d965bbSLois Curfman McInnes Input Parameter: 90304d965bbSLois Curfman McInnes . snes - the SNES context 90404d965bbSLois Curfman McInnes 90504d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 90604d965bbSLois Curfman McInnes */ 9075615d1e5SSatish Balay #undef __FUNC__ 908b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 909f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 9105e42470aSBarry Smith { 9116831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 912d132466eSBarry Smith int ierr; 9136daaf66cSBarry Smith 9143a40ed3dSBarry Smith PetscFunctionBegin; 9156831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 916d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 917d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 918d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 919d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 9203a40ed3dSBarry Smith PetscFunctionReturn(0); 9215e42470aSBarry Smith } 92204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92304d965bbSLois Curfman McInnes /* 9246831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 92504d965bbSLois Curfman McInnes 92604d965bbSLois Curfman McInnes Input Parameters: 92704d965bbSLois Curfman McInnes . SNES - the SNES context 92804d965bbSLois Curfman McInnes . viewer - visualization context 92904d965bbSLois Curfman McInnes 93004d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 93104d965bbSLois Curfman McInnes */ 9325615d1e5SSatish Balay #undef __FUNC__ 933b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 934e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 935a935fc98SLois Curfman McInnes { 9366831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 93719bcc07fSBarry Smith char *cstr; 93851695f54SLois Curfman McInnes int ierr; 9396831982aSBarry Smith PetscTruth isascii; 940a935fc98SLois Curfman McInnes 9413a40ed3dSBarry Smith PetscFunctionBegin; 9426831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 9430f5bd95cSBarry Smith if (isascii) { 94419bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 94519bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 94619bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 94719bcc07fSBarry Smith else cstr = "unknown"; 948d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 949d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9505cd90555SBarry Smith } else { 9510f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 95219bcc07fSBarry Smith } 9533a40ed3dSBarry Smith PetscFunctionReturn(0); 954a935fc98SLois Curfman McInnes } 95504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 95604d965bbSLois Curfman McInnes /* 9576831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 95804d965bbSLois Curfman McInnes 95904d965bbSLois Curfman McInnes Input Parameter: 96004d965bbSLois Curfman McInnes . snes - the SNES context 96104d965bbSLois Curfman McInnes 96204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 96304d965bbSLois Curfman McInnes */ 9645615d1e5SSatish Balay #undef __FUNC__ 965b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 966f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 9675e42470aSBarry Smith { 9686831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 9695e42470aSBarry Smith char ver[16]; 9705e42470aSBarry Smith double tmp; 971f1af5d2fSBarry Smith int ierr; 972f1af5d2fSBarry Smith PetscTruth flg; 9735e42470aSBarry Smith 9743a40ed3dSBarry Smith PetscFunctionBegin; 97509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 97625cdf11fSBarry Smith if (flg) { 9775e42470aSBarry Smith ls->alpha = tmp; 9785e42470aSBarry Smith } 97909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 98025cdf11fSBarry Smith if (flg) { 9815e42470aSBarry Smith ls->maxstep = tmp; 9825e42470aSBarry Smith } 98309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 98425cdf11fSBarry Smith if (flg) { 9855e42470aSBarry Smith ls->steptol = tmp; 9865e42470aSBarry Smith } 98709d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 98825cdf11fSBarry Smith if (flg) { 989c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9900f5bd95cSBarry Smith 991c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 992c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 993c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 994c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9950f5bd95cSBarry Smith 9960f5bd95cSBarry Smith if (isbasic) { 997af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9980f5bd95cSBarry Smith } else if (isnonorms) { 999af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 10000f5bd95cSBarry Smith } else if (isquad) { 1001af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 10020f5bd95cSBarry Smith } else if (iscubic) { 1003af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 10045e42470aSBarry Smith } 1005a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 10065e42470aSBarry Smith } 10073a40ed3dSBarry Smith PetscFunctionReturn(0); 10085e42470aSBarry Smith } 100904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 101004d965bbSLois Curfman McInnes /* 10116831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 10126831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 101304d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 101404d965bbSLois Curfman McInnes 101504d965bbSLois Curfman McInnes 101604d965bbSLois Curfman McInnes Input Parameter: 101704d965bbSLois Curfman McInnes . snes - the SNES context 101804d965bbSLois Curfman McInnes 101904d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 102004d965bbSLois Curfman McInnes */ 1021fb2e594dSBarry Smith EXTERN_C_BEGIN 10225615d1e5SSatish Balay #undef __FUNC__ 1023b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 1024f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 10255e42470aSBarry Smith { 102682bf6240SBarry Smith int ierr; 10276831982aSBarry Smith SNES_EQ_LS *neP; 10285e42470aSBarry Smith 10293a40ed3dSBarry Smith PetscFunctionBegin; 1030a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1031a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1032a8c6a408SBarry Smith } 103382bf6240SBarry Smith 1034f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 1035f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 1036f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 1037f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 1038f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 1039f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 1040f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 10415baf8537SBarry Smith snes->nwork = 0; 10425e42470aSBarry Smith 10436831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 10446831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 10455e42470aSBarry Smith snes->data = (void*)neP; 10465e42470aSBarry Smith neP->alpha = 1.e-4; 10475e42470aSBarry Smith neP->maxstep = 1.e8; 10485e42470aSBarry Smith neP->steptol = 1.e-12; 10495e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1050c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1051c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1052c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 105382bf6240SBarry Smith 10545d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 10555d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 105682bf6240SBarry Smith 10573a40ed3dSBarry Smith PetscFunctionReturn(0); 10585e42470aSBarry Smith } 1059fb2e594dSBarry Smith EXTERN_C_END 106082bf6240SBarry Smith 106182bf6240SBarry Smith 106282bf6240SBarry Smith 1063