163dd3a1aSKris Buschelman #define PETSCSNES_DLL 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" 13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 16dfbe8321SBarry Smith PetscErrorCode 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); 2663ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm));CHKERRQ(ierr); 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 33abb0e124SMatthew Knepley ierr = VecSetRandom(W,PETSC_NULL);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); 4063ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1));CHKERRQ(ierr); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 475ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4874637425SBarry Smith */ 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5274637425SBarry Smith { 5374637425SBarry Smith PetscReal a1,a2; 54dfbe8321SBarry Smith PetscErrorCode ierr; 5574637425SBarry Smith PetscTruth hastranspose; 5674637425SBarry Smith 5774637425SBarry Smith PetscFunctionBegin; 5874637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 5974637425SBarry Smith if (hastranspose) { 6074637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6179f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6274637425SBarry Smith 6374637425SBarry Smith /* Compute || J^T W|| */ 6474637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6574637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 679e247f21SBarry Smith if (a1) { 6863ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr); 6985471664SBarry Smith } 7074637425SBarry Smith } 7174637425SBarry Smith PetscFunctionReturn(0); 7274637425SBarry Smith } 7374637425SBarry Smith 7404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7504d965bbSLois Curfman McInnes 7604d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7794b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7804d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7904d965bbSLois Curfman McInnes respectively. 8004d965bbSLois Curfman McInnes 81fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8204d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8304d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 84fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8504d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8604d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 874b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 884b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 8904d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9004d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9104d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9204d965bbSLois Curfman McInnes for all nonlinear solvers. 9304d965bbSLois Curfman McInnes 9404d965bbSLois Curfman McInnes Another key routine is: 9504d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9604d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9704d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9804d965bbSLois Curfman McInnes SNESSolve() if necessary. 9904d965bbSLois Curfman McInnes 10004d965bbSLois Curfman McInnes Additional basic routines are: 10104d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10204d965bbSLois Curfman McInnes have actually been used. 10304d965bbSLois Curfman McInnes These are called by application codes via the interface routines 104186905e3SBarry Smith SNESView(). 10504d965bbSLois Curfman McInnes 10604d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10704d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10804d965bbSLois Curfman McInnes above description applies to these categories also. 10904d965bbSLois Curfman McInnes 11004d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1115e42470aSBarry Smith /* 1124b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11304d965bbSLois Curfman McInnes method with a line search. 1145e76c431SBarry Smith 11504d965bbSLois Curfman McInnes Input Parameters: 11604d965bbSLois Curfman McInnes . snes - the SNES context 1175e76c431SBarry Smith 11804d965bbSLois Curfman McInnes Output Parameter: 119c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12004d965bbSLois Curfman McInnes 12104d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1225e76c431SBarry Smith 1235e76c431SBarry Smith Notes: 1245e76c431SBarry Smith This implements essentially a truncated Newton method with a 1255e76c431SBarry Smith line search. By default a cubic backtracking line search 1265e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1275e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 128393d2d9aSLois Curfman McInnes and Schnabel. 1295e76c431SBarry Smith */ 1304a2ae208SSatish Balay #undef __FUNCT__ 1314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 132dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1335e76c431SBarry Smith { 1344b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1356849ba73SBarry Smith PetscErrorCode ierr; 136a7cc72afSBarry Smith PetscInt maxits,i,lits; 137a7cc72afSBarry Smith PetscTruth lssucceed; 138112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1405e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 141c293cc10SBarry Smith KSP ksp; 1425e76c431SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 14494b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 145184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 146184914b5SBarry Smith 1475e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1485e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1505e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1515e42470aSBarry Smith G = snes->work[1]; 1525e42470aSBarry Smith W = snes->work[2]; 1535e76c431SBarry Smith 1544c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1554c49b128SBarry Smith snes->iter = 0; 1564c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15719717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 158cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 15943e71028SBarry Smith if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1600f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1615e42470aSBarry Smith snes->norm = fnorm; 1620f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16384c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16494a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1655e76c431SBarry Smith 16670441072SBarry Smith if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16794a9d846SBarry Smith 168d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 169d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 170d034289bSBarry Smith 1715e76c431SBarry Smith for (i=0; i<maxits; i++) { 1725e76c431SBarry Smith 173329e7e40SMatthew Knepley /* Call general purpose update function */ 174abc0a331SBarry Smith if (snes->update) { 175329e7e40SMatthew Knepley ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr); 176de8cb200SMatthew Knepley } 177329e7e40SMatthew Knepley 178ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1795334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 18094b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 18123ce1328SBarry Smith ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 182c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 18374637425SBarry Smith 1849c3ca13aSBarry Smith if (neP->precheckstep) { 1859c3ca13aSBarry Smith PetscTruth changed_y = PETSC_FALSE; 1869c3ca13aSBarry Smith ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 1879c3ca13aSBarry Smith } 1889c3ca13aSBarry Smith 189b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 19074637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 19185471664SBarry Smith } 19274637425SBarry Smith 19390cbd584SBarry Smith /* should check what happened to the linear solve? */ 1943505fcc1SBarry Smith snes->linear_its += lits; 19563ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr); 196ea4d3ed3SLois Curfman McInnes 197ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 198ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 199ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 200ea4d3ed3SLois Curfman McInnes */ 20181b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 202a7cc72afSBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 20363ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed));CHKERRQ(ierr); 204a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 2053c632250SBarry Smith TMP = X; X = W; snes->vec_sol_always = X; W = TMP; 206a3b891d8SBarry Smith fnorm = gnorm; 207*4a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 208a3b891d8SBarry Smith 2095ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2105ed2d596SBarry Smith snes->iter = i+1; 2115ed2d596SBarry Smith snes->norm = fnorm; 2125ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2135ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2145ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2155ed2d596SBarry Smith 216a7cc72afSBarry Smith if (!lssucceed) { 2178a5d9ee4SBarry Smith PetscTruth ismin; 21850ffb88aSMatthew Knepley 21950ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2203505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2218a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2228a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2233505fcc1SBarry Smith break; 2243505fcc1SBarry Smith } 22550ffb88aSMatthew Knepley } 2265e76c431SBarry Smith 2275e76c431SBarry Smith /* Test for convergence */ 22829e0b56fSBarry Smith if (snes->converged) { 22929e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2303505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2313505fcc1SBarry Smith if (snes->reason) { 23216c95cacSBarry Smith break; 23316c95cacSBarry Smith } 23416c95cacSBarry Smith } 23529e0b56fSBarry Smith } 23639e2f89bSBarry Smith if (X != snes->vec_sol) { 237393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 238e15f7bb5SBarry Smith } 239e15f7bb5SBarry Smith if (F != snes->vec_func) { 240e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 241e15f7bb5SBarry Smith } 24239e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 24339e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 24452392280SLois Curfman McInnes if (i == maxits) { 24563ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr); 2463505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 24752392280SLois Curfman McInnes } 2483a40ed3dSBarry Smith PetscFunctionReturn(0); 2495e76c431SBarry Smith } 25004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25104d965bbSLois Curfman McInnes /* 2524b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2534b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 25404d965bbSLois Curfman McInnes 25504d965bbSLois Curfman McInnes Input Parameter: 25604d965bbSLois Curfman McInnes . snes - the SNES context 25704d965bbSLois Curfman McInnes . x - the solution vector 25804d965bbSLois Curfman McInnes 25904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26004d965bbSLois Curfman McInnes 26104d965bbSLois Curfman McInnes Notes: 2624b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 26304d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 26404d965bbSLois Curfman McInnes the call to SNESSolve(). 26504d965bbSLois Curfman McInnes */ 2664a2ae208SSatish Balay #undef __FUNCT__ 2674b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 268dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2695e76c431SBarry Smith { 270dfbe8321SBarry Smith PetscErrorCode ierr; 2713a40ed3dSBarry Smith 2723a40ed3dSBarry Smith PetscFunctionBegin; 273e74804ceSBarry Smith if (!snes->work) { 27481b6cf68SLois Curfman McInnes snes->nwork = 4; 275d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 276efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 27781b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 278e74804ceSBarry Smith } 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 2805e76c431SBarry Smith } 28104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 28204d965bbSLois Curfman McInnes /* 2834b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2844b27c08aSLois Curfman McInnes with SNESCreate_LS(). 28504d965bbSLois Curfman McInnes 28604d965bbSLois Curfman McInnes Input Parameter: 28704d965bbSLois Curfman McInnes . snes - the SNES context 28804d965bbSLois Curfman McInnes 289de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 29004d965bbSLois Curfman McInnes */ 2914a2ae208SSatish Balay #undef __FUNCT__ 2924b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 293dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 2945e76c431SBarry Smith { 295dfbe8321SBarry Smith PetscErrorCode ierr; 2963a40ed3dSBarry Smith 2973a40ed3dSBarry Smith PetscFunctionBegin; 2985baf8537SBarry Smith if (snes->nwork) { 2994b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 30021c89e3eSBarry Smith } 3015c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 3035e76c431SBarry Smith } 30404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3054a2ae208SSatish Balay #undef __FUNCT__ 30617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 30782bf6240SBarry Smith 3084b828684SBarry Smith /*@C 30917bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3105e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3115e76c431SBarry Smith to serve as a template and is not recommended for general use. 3125e76c431SBarry Smith 313c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 314c7afd0dbSLois Curfman McInnes 3155e76c431SBarry Smith Input Parameters: 316c7afd0dbSLois Curfman McInnes + snes - nonlinear context 317af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3185e76c431SBarry Smith . x - current iterate 3195e76c431SBarry Smith . f - residual evaluated at x 3203c632250SBarry Smith . y - search direction 3215e76c431SBarry Smith . w - work vector 322c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3235e76c431SBarry Smith 324c4a48953SLois Curfman McInnes Output Parameters: 325c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3263c632250SBarry Smith . w - new iterate 3275e76c431SBarry Smith . gnorm - 2-norm of g 3285e76c431SBarry Smith . ynorm - 2-norm of search length 329a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 330fee21e36SBarry Smith 331c4a48953SLois Curfman McInnes Options Database Key: 33217bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 333c4a48953SLois Curfman McInnes 33436851e7fSLois Curfman McInnes Level: advanced 33536851e7fSLois Curfman McInnes 33628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 33728ae5a21SLois Curfman McInnes 33817bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 33917bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3405e76c431SBarry Smith @*/ 34117bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 3425e76c431SBarry Smith { 343dfbe8321SBarry Smith PetscErrorCode ierr; 3444b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3453c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3465334005bSBarry Smith 3473a40ed3dSBarry Smith PetscFunctionBegin; 348a7cc72afSBarry Smith *flag = PETSC_TRUE; 349d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 350a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 35179f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3523c632250SBarry Smith if (neP->postcheckstep) { 3533c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3548f99978dSLois Curfman McInnes } 3553c632250SBarry Smith if (changed_y) { 35679f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3573c632250SBarry Smith } 3583c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 359d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 36019717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 36119717074SBarry Smith } 362d5e45103SBarry Smith CHKERRQ(ierr); 363d5e45103SBarry Smith 364a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 36543e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 366d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3673a40ed3dSBarry Smith PetscFunctionReturn(0); 3685e76c431SBarry Smith } 36904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3704a2ae208SSatish Balay #undef __FUNCT__ 37117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 37282bf6240SBarry Smith 37329e0b56fSBarry Smith /*@C 37417bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 37529e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 37629e0b56fSBarry Smith even compute the norm of the function or search direction; this 37729e0b56fSBarry Smith is intended only when you know the full step is fine and are 37829e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 379c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 380c7afd0dbSLois Curfman McInnes 381c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 38229e0b56fSBarry Smith 38329e0b56fSBarry Smith Input Parameters: 384c7afd0dbSLois Curfman McInnes + snes - nonlinear context 385af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 38629e0b56fSBarry Smith . x - current iterate 38729e0b56fSBarry Smith . f - residual evaluated at x 3883c632250SBarry Smith . y - search direction 38929e0b56fSBarry Smith . w - work vector 390c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 39129e0b56fSBarry Smith 39229e0b56fSBarry Smith Output Parameters: 393c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3943c632250SBarry Smith . w - new iterate 39529e0b56fSBarry Smith . gnorm - not changed 39629e0b56fSBarry Smith . ynorm - not changed 397a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 398fee21e36SBarry Smith 39929e0b56fSBarry Smith Options Database Key: 40017bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 40129e0b56fSBarry Smith 4028cbba510SBarry Smith Notes: 40317bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 404ea56f5baSLois Curfman McInnes either the options 405ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 406ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 407329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 408329f5518SBarry Smith otherwise, the SNES solver will generate an error. 409329f5518SBarry Smith 410329f5518SBarry Smith During the final iteration this will not evaluate the function at 411329f5518SBarry Smith the solution point. This is to save a function evaluation while 412329f5518SBarry Smith using pseudo-timestepping. 4138cbba510SBarry Smith 414ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 415ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 416ea56f5baSLois Curfman McInnes correct, since they are not computed. 417ea56f5baSLois Curfman McInnes 418ea56f5baSLois Curfman McInnes Level: advanced 4198cbba510SBarry Smith 42029e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 42129e0b56fSBarry Smith 42217bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 42317bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 42429e0b56fSBarry Smith @*/ 42517bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 42629e0b56fSBarry Smith { 427dfbe8321SBarry Smith PetscErrorCode ierr; 4284b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4293c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 43029e0b56fSBarry Smith 4313a40ed3dSBarry Smith PetscFunctionBegin; 432a7cc72afSBarry Smith *flag = PETSC_TRUE; 433d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 43479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4353c632250SBarry Smith if (neP->postcheckstep) { 4363c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4373c632250SBarry Smith } 4383c632250SBarry Smith if (changed_y) { 43979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4408f99978dSLois Curfman McInnes } 441329f5518SBarry Smith 442329f5518SBarry Smith /* don't evaluate function the last time through */ 443329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4443c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 445d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 44619717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 44719717074SBarry Smith } 448d5e45103SBarry Smith CHKERRQ(ierr); 449329f5518SBarry Smith } 450d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4513a40ed3dSBarry Smith PetscFunctionReturn(0); 45229e0b56fSBarry Smith } 45304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4544a2ae208SSatish Balay #undef __FUNCT__ 45517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4564b828684SBarry Smith /*@C 45717bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4585e76c431SBarry Smith 459c7afd0dbSLois Curfman McInnes Collective on SNES 460c7afd0dbSLois Curfman McInnes 4615e76c431SBarry Smith Input Parameters: 462c7afd0dbSLois Curfman McInnes + snes - nonlinear context 463af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4645e76c431SBarry Smith . x - current iterate 4655e76c431SBarry Smith . f - residual evaluated at x 4663c632250SBarry Smith . y - search direction 4675e76c431SBarry Smith . w - work vector 468c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4695e76c431SBarry Smith 470393d2d9aSLois Curfman McInnes Output Parameters: 471c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4723c632250SBarry Smith . w - new iterate 4735e76c431SBarry Smith . gnorm - 2-norm of g 4745e76c431SBarry Smith . ynorm - 2-norm of search length 475a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 476fee21e36SBarry Smith 477c4a48953SLois Curfman McInnes Options Database Key: 47817bae607SBarry Smith $ -snes_ls cubic - Activates SNESLineSearchCubic() 479c4a48953SLois Curfman McInnes 4805e76c431SBarry Smith Notes: 4815e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4825e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4835e76c431SBarry Smith 48436851e7fSLois Curfman McInnes Level: advanced 48536851e7fSLois Curfman McInnes 48628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 48728ae5a21SLois Curfman McInnes 48817bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 4895e76c431SBarry Smith @*/ 49017bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 4915e76c431SBarry Smith { 492406556e6SLois Curfman McInnes /* 493406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 494406556e6SLois Curfman McInnes minimization problem: 495406556e6SLois Curfman McInnes min z(x): R^n -> R, 496406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 497406556e6SLois Curfman McInnes */ 498406556e6SLois Curfman McInnes 4995ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 500329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 501aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 502ea709b57SSatish Balay PetscScalar cinitslope,clambda; 5036b5873e3SBarry Smith #endif 504dfbe8321SBarry Smith PetscErrorCode ierr; 505a7cc72afSBarry Smith PetscInt count; 5064b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 50779f36460SBarry Smith PetscScalar scale; 5083c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 5095e76c431SBarry Smith 5103a40ed3dSBarry Smith PetscFunctionBegin; 511d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 512a7cc72afSBarry Smith *flag = PETSC_TRUE; 5135e76c431SBarry Smith alpha = neP->alpha; 5145e76c431SBarry Smith maxstep = neP->maxstep; 5155e76c431SBarry Smith steptol = neP->steptol; 5165e76c431SBarry Smith 517cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 5189e247f21SBarry Smith if (!*ynorm) { 51917bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Search direction and size is 0\n"));CHKERRQ(ierr); 520a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5213c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 522a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 523a7cc72afSBarry Smith *flag = PETSC_FALSE; 524ad922ac9SBarry Smith goto theend1; 52594a9d846SBarry Smith } 5265e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5275e42470aSBarry Smith scale = maxstep/(*ynorm); 528aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 52917bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr); 5306b5873e3SBarry Smith #else 53117bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr); 5326b5873e3SBarry Smith #endif 5332dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 5345e76c431SBarry Smith *ynorm = maxstep; 5355e76c431SBarry Smith } 536ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5375ca10a37SBarry Smith minlambda = steptol/rellength; 538a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 539aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 540a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 541329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5425e42470aSBarry Smith #else 543a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5445e42470aSBarry Smith #endif 5455e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5465e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5475e76c431SBarry Smith 548393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 54979f36460SBarry Smith ierr = VecAYPX(w,-1.0,x);CHKERRQ(ierr); 55043e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 55117bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 55243e71028SBarry Smith *flag = PETSC_FALSE; 55343e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 55443e71028SBarry Smith goto theend1; 55543e71028SBarry Smith } 55619717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 557d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 55819717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 55919717074SBarry Smith } 560d5e45103SBarry Smith CHKERRQ(ierr); 561313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 56243e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5635d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 56417bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Using full step\n"));CHKERRQ(ierr); 565ad922ac9SBarry Smith goto theend1; 5665e76c431SBarry Smith } 5675e76c431SBarry Smith 5685e76c431SBarry Smith /* Fit points with quadratic */ 569313b5538SBarry Smith lambda = 1.0; 570a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5715e76c431SBarry Smith lambdaprev = lambda; 5725e76c431SBarry Smith gnormprev = *gnorm; 573329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 574ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 575ddd12767SBarry 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) 5792dcb1b2aSMatthew Knepley clambda = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr); 5805e42470aSBarry Smith #else 5812dcb1b2aSMatthew Knepley ierr = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr); 5825e42470aSBarry Smith #endif 58343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 58417bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr); 58543e71028SBarry Smith *flag = PETSC_FALSE; 58643e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 58743e71028SBarry Smith goto theend1; 58843e71028SBarry Smith } 58919717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 590d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 59119717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 59219717074SBarry Smith } 593d5e45103SBarry Smith CHKERRQ(ierr); 594cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 59543e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5965ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 59717bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 598ad922ac9SBarry Smith goto theend1; 5995e76c431SBarry Smith } 6005e76c431SBarry Smith 6015e76c431SBarry Smith /* Fit points with cubic */ 6025e76c431SBarry Smith count = 1; 6038229bfc2SMatthew Knepley while (count < 10000) { 6045e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 60517bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 60617bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 60743e71028SBarry Smith *flag = PETSC_FALSE; 60843e71028SBarry Smith break; 6095e76c431SBarry Smith } 6105d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6115d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 612ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6132b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6145e76c431SBarry Smith d = b*b - 3*a*initslope; 6155e76c431SBarry Smith if (d < 0.0) d = 0.0; 6165e76c431SBarry Smith if (a == 0.0) { 6175e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6185e76c431SBarry Smith } else { 6195e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6205e76c431SBarry Smith } 6215e76c431SBarry Smith lambdaprev = lambda; 6225e76c431SBarry Smith gnormprev = *gnorm; 623329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 624329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6255e76c431SBarry Smith else lambda = lambdatemp; 626393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 627ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 628aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 629ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 6302dcb1b2aSMatthew Knepley ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr); 6315e42470aSBarry Smith #else 6322dcb1b2aSMatthew Knepley ierr = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr); 6335e42470aSBarry Smith #endif 63443e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 63517bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 63617bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 637ed98166eSMatthew Knepley *flag = PETSC_FALSE; 63843e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 639ed98166eSMatthew Knepley break; 640ed98166eSMatthew Knepley } 64119717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 642d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 64319717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 64419717074SBarry Smith } 645d5e45103SBarry Smith CHKERRQ(ierr); 646cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 64743e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6485ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 64917bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 65043e71028SBarry Smith break; 6512b022350SLois Curfman McInnes } else { 65217bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 6535e76c431SBarry Smith } 6545e76c431SBarry Smith count++; 6555e76c431SBarry Smith } 6568229bfc2SMatthew Knepley if (count >= 10000) { 657cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6588229bfc2SMatthew Knepley } 659ad922ac9SBarry Smith theend1: 6608f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6613c632250SBarry Smith if (neP->postcheckstep && *flag) { 6623c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6633c632250SBarry Smith if (changed_y) { 66479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6653c632250SBarry Smith } 6663c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 6673c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 668d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 66919717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 67019717074SBarry Smith } 671d5e45103SBarry Smith CHKERRQ(ierr); 6728f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 67343e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6743c632250SBarry Smith ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 6758f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6763c632250SBarry Smith ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 6778f99978dSLois Curfman McInnes } 6788f99978dSLois Curfman McInnes } 679d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6803a40ed3dSBarry Smith PetscFunctionReturn(0); 6815e76c431SBarry Smith } 68204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6834a2ae208SSatish Balay #undef __FUNCT__ 68417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 6854b828684SBarry Smith /*@C 68617bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 6875e76c431SBarry Smith 688c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 689c7afd0dbSLois Curfman McInnes 6905e42470aSBarry Smith Input Parameters: 691c7afd0dbSLois Curfman McInnes + snes - the SNES context 692af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6935e42470aSBarry Smith . x - current iterate 6945e42470aSBarry Smith . f - residual evaluated at x 6953c632250SBarry Smith . y - search direction 6965e42470aSBarry Smith . w - work vector 697c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6985e42470aSBarry Smith 699c4a48953SLois Curfman McInnes Output Parameters: 700c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7013c632250SBarry Smith . w - new iterate 7025e42470aSBarry Smith . gnorm - 2-norm of g 7035e42470aSBarry Smith . ynorm - 2-norm of search length 704a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 705fee21e36SBarry Smith 706c4a48953SLois Curfman McInnes Options Database Key: 70717bae607SBarry Smith . -snes_ls quadratic - Activates SNESLineSearchQuadratic() 708c4a48953SLois Curfman McInnes 7095e42470aSBarry Smith Notes: 71017bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7115e42470aSBarry Smith 71236851e7fSLois Curfman McInnes Level: advanced 71336851e7fSLois Curfman McInnes 714f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 715f59ffdedSLois Curfman McInnes 71617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7175e42470aSBarry Smith @*/ 71817bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 7195e76c431SBarry Smith { 720406556e6SLois Curfman McInnes /* 721406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 722406556e6SLois Curfman McInnes minimization problem: 723406556e6SLois Curfman McInnes min z(x): R^n -> R, 724406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 725406556e6SLois Curfman McInnes */ 7265ca10a37SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 727aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 728ea709b57SSatish Balay PetscScalar cinitslope,clambda; 7296b5873e3SBarry Smith #endif 730dfbe8321SBarry Smith PetscErrorCode ierr; 731a7cc72afSBarry Smith PetscInt count; 7324b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 73379f36460SBarry Smith PetscScalar scale; 7343c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7355e76c431SBarry Smith 7363a40ed3dSBarry Smith PetscFunctionBegin; 737d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 738a7cc72afSBarry Smith *flag = PETSC_TRUE; 7395e42470aSBarry Smith alpha = neP->alpha; 7405e42470aSBarry Smith maxstep = neP->maxstep; 7415e42470aSBarry Smith steptol = neP->steptol; 7425e76c431SBarry Smith 7433505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 74463b9fa5eSBarry Smith if (*ynorm == 0.0) { 74517bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic: Search direction and size is 0\n"));CHKERRQ(ierr); 746b37302c6SLois Curfman McInnes *gnorm = fnorm; 747b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 748b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 749a7cc72afSBarry Smith *flag = PETSC_FALSE; 750ad922ac9SBarry Smith goto theend2; 75194a9d846SBarry Smith } 7525e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 7535e42470aSBarry Smith scale = maxstep/(*ynorm); 7542dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 7555e42470aSBarry Smith *ynorm = maxstep; 7565e76c431SBarry Smith } 757ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 7585ca10a37SBarry Smith minlambda = steptol/rellength; 759a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 760aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 761a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 762329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7635e42470aSBarry Smith #else 764a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7655e42470aSBarry Smith #endif 7665e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7675e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7685e42470aSBarry Smith 769393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 77079f36460SBarry Smith ierr = VecAYPX(w,-1.0,x);CHKERRQ(ierr); 77143e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 77217bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 77343e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 77443e71028SBarry Smith *flag = PETSC_FALSE; 77543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 77643e71028SBarry Smith goto theend2; 77743e71028SBarry Smith } 77819717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 779d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 78019717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 78119717074SBarry Smith } 782d5e45103SBarry Smith CHKERRQ(ierr); 783cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 78443e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 7855d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 786393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 78717bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic: Using full step\n"));CHKERRQ(ierr); 788ad922ac9SBarry Smith goto theend2; 7895e42470aSBarry Smith } 7905e42470aSBarry Smith 7915e42470aSBarry Smith /* Fit points with quadratic */ 792313b5538SBarry Smith lambda = 1.0; 7935e42470aSBarry Smith count = 1; 7945ca10a37SBarry Smith while (PETSC_TRUE) { 7955e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 79617bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 79717bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 798f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 79943e71028SBarry Smith *flag = PETSC_FALSE; 80043e71028SBarry Smith break; 8015e42470aSBarry Smith } 802a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 803329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 804329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 805329f5518SBarry Smith else lambda = lambdatemp; 806393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 8073505fcc1SBarry Smith lambdaneg = -lambda; 808aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 8092dcb1b2aSMatthew Knepley clambda = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr); 8105e42470aSBarry Smith #else 8112dcb1b2aSMatthew Knepley ierr = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr); 8125e42470aSBarry Smith #endif 81343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 81417bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 81517bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 81643e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 817ed98166eSMatthew Knepley *flag = PETSC_FALSE; 81843e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 819ed98166eSMatthew Knepley break; 820ed98166eSMatthew Knepley } 82119717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 822d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 82319717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 82419717074SBarry Smith } 825d5e45103SBarry Smith CHKERRQ(ierr); 826cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 82743e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8285ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 829393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 83017bae607SBarry Smith ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr); 83106259719SBarry Smith break; 8325e42470aSBarry Smith } 8335e42470aSBarry Smith count++; 8345e42470aSBarry Smith } 835ad922ac9SBarry Smith theend2: 8368f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8373c632250SBarry Smith if (neP->postcheckstep) { 8383c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8393c632250SBarry Smith if (changed_y) { 84079f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8413c632250SBarry Smith } 8423c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8433c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 844d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 84519717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 84619717074SBarry Smith } 847d5e45103SBarry Smith CHKERRQ(ierr); 8488f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 8493c632250SBarry Smith ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 8508f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 8513c632250SBarry Smith ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 8523c632250SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8538f99978dSLois Curfman McInnes } 8548f99978dSLois Curfman McInnes } 855d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8563a40ed3dSBarry Smith PetscFunctionReturn(0); 8575e76c431SBarry Smith } 8582343ba6eSBarry Smith 85904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8604a2ae208SSatish Balay #undef __FUNCT__ 86117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 862c9e6a524SLois Curfman McInnes /*@C 86317bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8644b27c08aSLois Curfman McInnes by the method SNESLS. 8655e76c431SBarry Smith 8665e76c431SBarry Smith Input Parameters: 867c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 868af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 869c7afd0dbSLois Curfman McInnes - func - pointer to int function 8705e76c431SBarry Smith 871fee21e36SBarry Smith Collective on SNES 872fee21e36SBarry Smith 873c4a48953SLois Curfman McInnes Available Routines: 87417bae607SBarry Smith + SNESLineSearchCubic() - default line search 87517bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 87617bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 87717bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8785e76c431SBarry Smith 879c4a48953SLois Curfman McInnes Options Database Keys: 8804b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8814b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 8824b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 8834b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 8843304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 8853304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 886c4a48953SLois Curfman McInnes 8875e76c431SBarry Smith Calling sequence of func: 888bd208895SLois Curfman McInnes .vb 889af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 890a7cc72afSBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 891bd208895SLois Curfman McInnes .ve 8925e76c431SBarry Smith 8935e76c431SBarry Smith Input parameters for func: 894c7afd0dbSLois Curfman McInnes + snes - nonlinear context 895af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8965e76c431SBarry Smith . x - current iterate 8975e76c431SBarry Smith . f - residual evaluated at x 8983c632250SBarry Smith . y - search direction 8995e76c431SBarry Smith . w - work vector 900c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9015e76c431SBarry Smith 9025e76c431SBarry Smith Output parameters for func: 903c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9043c632250SBarry Smith . w - new iterate 9055e76c431SBarry Smith . gnorm - 2-norm of g 9065e76c431SBarry Smith . ynorm - 2-norm of search length 907a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 908f59ffdedSLois Curfman McInnes 90936851e7fSLois Curfman McInnes Level: advanced 91036851e7fSLois Curfman McInnes 911f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 912f59ffdedSLois Curfman McInnes 91317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 91417bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9155e76c431SBarry Smith @*/ 91617bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 9175e76c431SBarry Smith { 918a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 91982bf6240SBarry Smith 9203a40ed3dSBarry Smith PetscFunctionBegin; 92117bae607SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 92282bf6240SBarry Smith if (f) { 923af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 92482bf6240SBarry Smith } 9253a40ed3dSBarry Smith PetscFunctionReturn(0); 9265e76c431SBarry Smith } 9278e019c35SBarry Smith 928a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 92904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 930fb2e594dSBarry Smith EXTERN_C_BEGIN 9314a2ae208SSatish Balay #undef __FUNCT__ 93217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 93317bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 93482bf6240SBarry Smith { 93582bf6240SBarry Smith PetscFunctionBegin; 9364b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9374b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 93882bf6240SBarry Smith PetscFunctionReturn(0); 93982bf6240SBarry Smith } 940fb2e594dSBarry Smith EXTERN_C_END 94104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9424a2ae208SSatish Balay #undef __FUNCT__ 9433c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 944c8dd0c0dSLois Curfman McInnes /*@C 9453c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9464b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 947c8dd0c0dSLois Curfman McInnes 948c8dd0c0dSLois Curfman McInnes Input Parameters: 949c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9503c632250SBarry Smith . func - pointer to function 951c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 952c8dd0c0dSLois Curfman McInnes 953c8dd0c0dSLois Curfman McInnes Collective on SNES 954c8dd0c0dSLois Curfman McInnes 955c8dd0c0dSLois Curfman McInnes Calling sequence of func: 956c8dd0c0dSLois Curfman McInnes .vb 9573c632250SBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 958c8dd0c0dSLois Curfman McInnes .ve 959b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 960b1ae27eaSLois Curfman McInnes on failure. 961c8dd0c0dSLois Curfman McInnes 962c8dd0c0dSLois Curfman McInnes Input parameters for func: 963c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 964c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9653c632250SBarry Smith . x - previous iterate 9663c632250SBarry Smith . y - new search direction and length 9673c632250SBarry Smith - w - current candidate iterate 968c8dd0c0dSLois Curfman McInnes 969c8dd0c0dSLois Curfman McInnes Output parameters for func: 9703c632250SBarry Smith + y - search direction (possibly changed) 9713c632250SBarry Smith . w - current iterate (possibly modified) 9723c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 9733c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 974c8dd0c0dSLois Curfman McInnes 975c8dd0c0dSLois Curfman McInnes Level: advanced 976c8dd0c0dSLois Curfman McInnes 9779e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 9789e247f21SBarry Smith 9793c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 9803c632250SBarry Smith 9813c632250SBarry Smith On input w = x + y 9823c632250SBarry Smith 98317bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 984b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 985ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9868f99978dSLois Curfman McInnes 98717bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 988ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 989ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 990ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 9919e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 992b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9938f99978dSLois Curfman McInnes 994c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 995c8dd0c0dSLois Curfman McInnes 99617bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 997c8dd0c0dSLois Curfman McInnes @*/ 9983c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 999c8dd0c0dSLois Curfman McInnes { 10003c632250SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 1001c8dd0c0dSLois Curfman McInnes 1002c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10033c632250SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1004c8dd0c0dSLois Curfman McInnes if (f) { 1005c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1006c8dd0c0dSLois Curfman McInnes } 1007c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1008c8dd0c0dSLois Curfman McInnes } 10099c3ca13aSBarry Smith #undef __FUNCT__ 10109c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10119c3ca13aSBarry Smith /*@C 10129c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10139c3ca13aSBarry Smith 10149c3ca13aSBarry Smith Input Parameters: 10159c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10169c3ca13aSBarry Smith . func - pointer to function 10179c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10189c3ca13aSBarry Smith 10199c3ca13aSBarry Smith Collective on SNES 10209c3ca13aSBarry Smith 10219c3ca13aSBarry Smith Calling sequence of func: 10229c3ca13aSBarry Smith .vb 10239c3ca13aSBarry Smith int func (SNES snes, Vec x,Vec y,,void *checkctx, PetscTruth *changed_y) 10249c3ca13aSBarry Smith .ve 10259c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10269c3ca13aSBarry Smith on failure. 10279c3ca13aSBarry Smith 10289c3ca13aSBarry Smith Input parameters for func: 10299c3ca13aSBarry Smith + snes - nonlinear context 10309c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10319c3ca13aSBarry Smith . x - previous iterate 10329c3ca13aSBarry Smith - y - new search direction and length 10339c3ca13aSBarry Smith 10349c3ca13aSBarry Smith Output parameters for func: 10359c3ca13aSBarry Smith + y - search direction (possibly changed) 10369c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10379c3ca13aSBarry Smith 10389c3ca13aSBarry Smith Level: advanced 10399c3ca13aSBarry Smith 10409c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10419c3ca13aSBarry Smith 10429c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10439c3ca13aSBarry Smith 104417bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck() 10459c3ca13aSBarry Smith @*/ 10469c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 10479c3ca13aSBarry Smith { 10489c3ca13aSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 10499c3ca13aSBarry Smith 10509c3ca13aSBarry Smith PetscFunctionBegin; 10519c3ca13aSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 10529c3ca13aSBarry Smith if (f) { 10539c3ca13aSBarry Smith ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 10549c3ca13aSBarry Smith } 10559c3ca13aSBarry Smith PetscFunctionReturn(0); 10569c3ca13aSBarry Smith } 10579c3ca13aSBarry Smith 1058c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 10599c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 10609c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1061c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 10624a2ae208SSatish Balay #undef __FUNCT__ 10633c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 10649c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1065c8dd0c0dSLois Curfman McInnes { 1066c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10673c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 10683c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1069c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1070c8dd0c0dSLois Curfman McInnes } 1071c8dd0c0dSLois Curfman McInnes EXTERN_C_END 10729c3ca13aSBarry Smith 10739c3ca13aSBarry Smith EXTERN_C_BEGIN 10749c3ca13aSBarry Smith #undef __FUNCT__ 10759c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 10769c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 10779c3ca13aSBarry Smith { 10789c3ca13aSBarry Smith PetscFunctionBegin; 10799c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 10809c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 10819c3ca13aSBarry Smith PetscFunctionReturn(0); 10829c3ca13aSBarry Smith } 10839c3ca13aSBarry Smith EXTERN_C_END 1084c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 108504d965bbSLois Curfman McInnes /* 10864b27c08aSLois Curfman McInnes SNESPrintHelp_LS - Prints all options for the SNES_LS method. 1087329e7e40SMatthew Knepley 1088329e7e40SMatthew Knepley Input Parameter: 1089329e7e40SMatthew Knepley . snes - the SNES context 1090329e7e40SMatthew Knepley 1091329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 1092329e7e40SMatthew Knepley */ 1093329e7e40SMatthew Knepley #undef __FUNCT__ 10944b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS" 10956849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p) 1096329e7e40SMatthew Knepley { 10974b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1098329e7e40SMatthew Knepley 1099329e7e40SMatthew Knepley PetscFunctionBegin; 11004b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 11014b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 11024b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 11034b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 11044b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 1105329e7e40SMatthew Knepley PetscFunctionReturn(0); 1106329e7e40SMatthew Knepley } 1107329e7e40SMatthew Knepley 1108329e7e40SMatthew Knepley /* 11094b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 111004d965bbSLois Curfman McInnes 111104d965bbSLois Curfman McInnes Input Parameters: 111204d965bbSLois Curfman McInnes . SNES - the SNES context 111304d965bbSLois Curfman McInnes . viewer - visualization context 111404d965bbSLois Curfman McInnes 111504d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 111604d965bbSLois Curfman McInnes */ 11174a2ae208SSatish Balay #undef __FUNCT__ 11184b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11196849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1120a935fc98SLois Curfman McInnes { 11214b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11222fc52814SBarry Smith const char *cstr; 1123dfbe8321SBarry Smith PetscErrorCode ierr; 112432077d6dSBarry Smith PetscTruth iascii; 1125a935fc98SLois Curfman McInnes 11263a40ed3dSBarry Smith PetscFunctionBegin; 112732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 112832077d6dSBarry Smith if (iascii) { 112917bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 113017bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 113117bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 113219bcc07fSBarry Smith else cstr = "unknown"; 1133b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1134b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 11355cd90555SBarry Smith } else { 113679a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 113719bcc07fSBarry Smith } 11383a40ed3dSBarry Smith PetscFunctionReturn(0); 1139a935fc98SLois Curfman McInnes } 114004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 114104d965bbSLois Curfman McInnes /* 11424b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 114304d965bbSLois Curfman McInnes 114404d965bbSLois Curfman McInnes Input Parameter: 114504d965bbSLois Curfman McInnes . snes - the SNES context 114604d965bbSLois Curfman McInnes 114704d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 114804d965bbSLois Curfman McInnes */ 11494a2ae208SSatish Balay #undef __FUNCT__ 11504b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11516849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11525e42470aSBarry Smith { 11534b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1154e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1155dfbe8321SBarry Smith PetscErrorCode ierr; 1156a7cc72afSBarry Smith PetscInt indx; 1157f1af5d2fSBarry Smith PetscTruth flg; 11585e42470aSBarry Smith 11593a40ed3dSBarry Smith PetscFunctionBegin; 1160b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 11614b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 11624b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 11634b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1164186905e3SBarry Smith 116517bae607SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 116625cdf11fSBarry Smith if (flg) { 116722e36657SBarry Smith switch (indx) { 1168b49fd9e1SBarry Smith case 0: 116917bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1170b49fd9e1SBarry Smith break; 1171b49fd9e1SBarry Smith case 1: 117217bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1173b49fd9e1SBarry Smith break; 1174b49fd9e1SBarry Smith case 2: 117517bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1176b49fd9e1SBarry Smith break; 1177b49fd9e1SBarry Smith case 3: 117817bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1179b49fd9e1SBarry Smith break; 11805e42470aSBarry Smith } 11815e42470aSBarry Smith } 1182b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11833a40ed3dSBarry Smith PetscFunctionReturn(0); 11845e42470aSBarry Smith } 118504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1186ebe3b25bSBarry Smith /*MC 1187ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 118804d965bbSLois Curfman McInnes 1189ebe3b25bSBarry Smith Options Database: 1190ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1191ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1192ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1193ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1194ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1195ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 119604d965bbSLois Curfman McInnes 1197ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 119804d965bbSLois Curfman McInnes 1199ee3001cbSBarry Smith Level: beginner 1200ee3001cbSBarry Smith 120117bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 120217bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 120317bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1204ebe3b25bSBarry Smith 1205ebe3b25bSBarry Smith M*/ 1206fb2e594dSBarry Smith EXTERN_C_BEGIN 12074a2ae208SSatish Balay #undef __FUNCT__ 12084b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 120963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 12105e42470aSBarry Smith { 1211dfbe8321SBarry Smith PetscErrorCode ierr; 12124b27c08aSLois Curfman McInnes SNES_LS *neP; 12135e42470aSBarry Smith 12143a40ed3dSBarry Smith PetscFunctionBegin; 12154b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 12164b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 12174b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 12184b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 12194b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 12204b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 12214b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 12225baf8537SBarry Smith snes->nwork = 0; 12235e42470aSBarry Smith 12244b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 122552e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 12265e42470aSBarry Smith snes->data = (void*)neP; 12275e42470aSBarry Smith neP->alpha = 1.e-4; 12285e42470aSBarry Smith neP->maxstep = 1.e8; 12295e42470aSBarry Smith neP->steptol = 1.e-12; 123017bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1231c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12323c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12333c632250SBarry Smith neP->postcheck = PETSC_NULL; 12343c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12353c632250SBarry Smith neP->precheck = PETSC_NULL; 123682bf6240SBarry Smith 123717bae607SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 12383c632250SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 12399c3ca13aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 124082bf6240SBarry Smith 12413a40ed3dSBarry Smith PetscFunctionReturn(0); 12425e42470aSBarry Smith } 1243fb2e594dSBarry Smith EXTERN_C_END 124482bf6240SBarry Smith 124582bf6240SBarry Smith 124682bf6240SBarry Smith 1247