163dd3a1aSKris Buschelman #define PETSCSNES_DLL 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 58a5d9ee4SBarry Smith /* 620f52e01SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 720f52e01SBarry Smith || F(u) ||_2 but not a zero, F(u) = 0. 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 920f52e01SBarry Smith for this trick. One assumes that the probability that W is in the null space of J is very, very small. 108a5d9ee4SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 131e2582c4SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,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); 261e2582c4SBarry Smith ierr = PetscInfo1(snes,"|| 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); 401e2582c4SBarry Smith ierr = PetscInfo1(snes,"(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" 511e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,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) { 681e2582c4SBarry Smith ierr = PetscInfo1(snes,"||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; 13985385478SLisandro Dalcin PetscReal fnorm,gnorm,xnorm=0,ynorm; 14085385478SLisandro Dalcin Vec Y,X,F,G,W; 1413d4c4710SBarry Smith KSPConvergedReason kspreason; 1425e76c431SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 1443d4c4710SBarry Smith snes->numFailures = 0; 1453d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 146184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 147184914b5SBarry Smith 1485e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1495e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15039e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1515e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1525e42470aSBarry Smith G = snes->work[1]; 1535e42470aSBarry Smith W = snes->work[2]; 1545e76c431SBarry Smith 1554c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1564c49b128SBarry Smith snes->iter = 0; 15785385478SLisandro Dalcin snes->norm = 0; 1584c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15919717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1604936397dSBarry Smith if (snes->domainerror) { 1614936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1624936397dSBarry Smith PetscFunctionReturn(0); 1634936397dSBarry Smith } 1642613ca53SBarry Smith ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1652613ca53SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1662613ca53SBarry Smith ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 1672613ca53SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 168f66fdb6dSSatish Balay if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1690f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1705e42470aSBarry Smith snes->norm = fnorm; 1710f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 17284c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 17394a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1743f149594SLisandro Dalcin 1753f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1763f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 17785385478SLisandro Dalcin /* test convergence */ 17885385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 17906ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 180d034289bSBarry Smith 1815e76c431SBarry Smith for (i=0; i<maxits; i++) { 1825e76c431SBarry Smith 183329e7e40SMatthew Knepley /* Call general purpose update function */ 184e7788613SBarry Smith if (snes->ops->update) { 185e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 186de8cb200SMatthew Knepley } 187329e7e40SMatthew Knepley 188ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1895334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 19094b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 19171f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 1923d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1933d4c4710SBarry Smith if (kspreason < 0) { 1943d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 195ef998cc9SBarry Smith ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 1963d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 197ab81109fSSatish Balay break; 1983d4c4710SBarry Smith } 1993d4c4710SBarry Smith } 2003d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2013f149594SLisandro Dalcin snes->linear_its += lits; 2023f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 20374637425SBarry Smith 2049c3ca13aSBarry Smith if (neP->precheckstep) { 2059c3ca13aSBarry Smith PetscTruth changed_y = PETSC_FALSE; 2069c3ca13aSBarry Smith ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 2079c3ca13aSBarry Smith } 2089c3ca13aSBarry Smith 209b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 2101e2582c4SBarry Smith ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 21185471664SBarry Smith } 21274637425SBarry Smith 213ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 214ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 215e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 216ea4d3ed3SLois Curfman McInnes */ 21785385478SLisandro Dalcin ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 2183f149594SLisandro Dalcin ynorm = 1; gnorm = fnorm; 219dc357ad2SBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 220ae15b995SBarry Smith ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 2214a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2224936397dSBarry Smith if (snes->domainerror) { 2234936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2244936397dSBarry Smith PetscFunctionReturn(0); 2254936397dSBarry Smith } 226a7cc72afSBarry Smith if (!lssucceed) { 22750ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 22885385478SLisandro Dalcin PetscTruth ismin; 2293505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2301e2582c4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 2318a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2323505fcc1SBarry Smith break; 2333505fcc1SBarry Smith } 23450ffb88aSMatthew Knepley } 23585385478SLisandro Dalcin /* Update function and solution vectors */ 23685385478SLisandro Dalcin fnorm = gnorm; 23785385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 23885385478SLisandro Dalcin ierr = VecCopy(W,X);CHKERRQ(ierr); 23985385478SLisandro Dalcin /* Monitor convergence */ 24085385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 24185385478SLisandro Dalcin snes->iter = i+1; 24285385478SLisandro Dalcin snes->norm = fnorm; 24385385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24485385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 24585385478SLisandro Dalcin SNESMonitor(snes,snes->iter,snes->norm); 24685385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 24785385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 248e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2493f149594SLisandro Dalcin if (snes->reason) break; 25016c95cacSBarry Smith } 25152392280SLois Curfman McInnes if (i == maxits) { 252ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25385385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25452392280SLois Curfman McInnes } 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 2565e76c431SBarry Smith } 25704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25804d965bbSLois Curfman McInnes /* 2594b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2604b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26104d965bbSLois Curfman McInnes 26204d965bbSLois Curfman McInnes Input Parameter: 26304d965bbSLois Curfman McInnes . snes - the SNES context 26404d965bbSLois Curfman McInnes . x - the solution vector 26504d965bbSLois Curfman McInnes 26604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26704d965bbSLois Curfman McInnes 26804d965bbSLois Curfman McInnes Notes: 2694b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 27004d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27104d965bbSLois Curfman McInnes the call to SNESSolve(). 27204d965bbSLois Curfman McInnes */ 2734a2ae208SSatish Balay #undef __FUNCT__ 2744b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 275dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2765e76c431SBarry Smith { 277dfbe8321SBarry Smith PetscErrorCode ierr; 2783a40ed3dSBarry Smith 2793a40ed3dSBarry Smith PetscFunctionBegin; 28085385478SLisandro Dalcin if (!snes->vec_sol_update) { 28185385478SLisandro Dalcin ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 28285385478SLisandro Dalcin ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 28385385478SLisandro Dalcin } 284e74804ceSBarry Smith if (!snes->work) { 28585385478SLisandro Dalcin snes->nwork = 3; 28685385478SLisandro Dalcin ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 287efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 288e74804ceSBarry Smith } 2893a40ed3dSBarry Smith PetscFunctionReturn(0); 2905e76c431SBarry Smith } 29104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 29204d965bbSLois Curfman McInnes /* 2934b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2944b27c08aSLois Curfman McInnes with SNESCreate_LS(). 29504d965bbSLois Curfman McInnes 29604d965bbSLois Curfman McInnes Input Parameter: 29704d965bbSLois Curfman McInnes . snes - the SNES context 29804d965bbSLois Curfman McInnes 299de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30004d965bbSLois Curfman McInnes */ 3014a2ae208SSatish Balay #undef __FUNCT__ 3024b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 303dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3045e76c431SBarry Smith { 305dfbe8321SBarry Smith PetscErrorCode ierr; 3063a40ed3dSBarry Smith 3073a40ed3dSBarry Smith PetscFunctionBegin; 30885385478SLisandro Dalcin if (snes->vec_sol_update) { 30985385478SLisandro Dalcin ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 31085385478SLisandro Dalcin snes->vec_sol_update = PETSC_NULL; 31185385478SLisandro Dalcin } 3125baf8537SBarry Smith if (snes->nwork) { 3134b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 31485385478SLisandro Dalcin snes->nwork = 0; 31585385478SLisandro Dalcin snes->work = PETSC_NULL; 31621c89e3eSBarry Smith } 3175c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 318557d3f75SLisandro Dalcin 319557d3f75SLisandro Dalcin /* clear composed functions */ 320557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 321557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 322557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 323557d3f75SLisandro Dalcin 3243a40ed3dSBarry Smith PetscFunctionReturn(0); 3255e76c431SBarry Smith } 32604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3274a2ae208SSatish Balay #undef __FUNCT__ 32817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 32982bf6240SBarry Smith 3304b828684SBarry Smith /*@C 33117bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3325e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3335e76c431SBarry Smith to serve as a template and is not recommended for general use. 3345e76c431SBarry Smith 335c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 336c7afd0dbSLois Curfman McInnes 3375e76c431SBarry Smith Input Parameters: 338c7afd0dbSLois Curfman McInnes + snes - nonlinear context 339af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3405e76c431SBarry Smith . x - current iterate 3415e76c431SBarry Smith . f - residual evaluated at x 3423c632250SBarry Smith . y - search direction 3435e76c431SBarry Smith . w - work vector 344dc357ad2SBarry Smith . fnorm - 2-norm of f 345dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 3465e76c431SBarry Smith 347c4a48953SLois Curfman McInnes Output Parameters: 348c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3493c632250SBarry Smith . w - new iterate 3505e76c431SBarry Smith . gnorm - 2-norm of g 3515e76c431SBarry Smith . ynorm - 2-norm of search length 352a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 353fee21e36SBarry Smith 354c4a48953SLois Curfman McInnes Options Database Key: 35517bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 356c4a48953SLois Curfman McInnes 35736851e7fSLois Curfman McInnes Level: advanced 35836851e7fSLois Curfman McInnes 35928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 36028ae5a21SLois Curfman McInnes 36117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 36217bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3635e76c431SBarry Smith @*/ 364dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 3655e76c431SBarry Smith { 366dfbe8321SBarry Smith PetscErrorCode ierr; 3674b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3683c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3695334005bSBarry Smith 3703a40ed3dSBarry Smith PetscFunctionBegin; 371a7cc72afSBarry Smith *flag = PETSC_TRUE; 372d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 373a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 37479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3753c632250SBarry Smith if (neP->postcheckstep) { 3763c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3778f99978dSLois Curfman McInnes } 3783c632250SBarry Smith if (changed_y) { 37979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3803c632250SBarry Smith } 3814936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 3824936397dSBarry Smith if (!snes->domainerror) { 383a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 384f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3854936397dSBarry Smith } 386d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3885e76c431SBarry Smith } 38904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3904a2ae208SSatish Balay #undef __FUNCT__ 39117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 39282bf6240SBarry Smith 39329e0b56fSBarry Smith /*@C 39417bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 39529e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 39629e0b56fSBarry Smith even compute the norm of the function or search direction; this 39729e0b56fSBarry Smith is intended only when you know the full step is fine and are 39829e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 399c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 400c7afd0dbSLois Curfman McInnes 401c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 40229e0b56fSBarry Smith 40329e0b56fSBarry Smith Input Parameters: 404c7afd0dbSLois Curfman McInnes + snes - nonlinear context 405af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 40629e0b56fSBarry Smith . x - current iterate 40729e0b56fSBarry Smith . f - residual evaluated at x 4083c632250SBarry Smith . y - search direction 40929e0b56fSBarry Smith . w - work vector 410dc357ad2SBarry Smith . fnorm - 2-norm of f 411dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 41229e0b56fSBarry Smith 41329e0b56fSBarry Smith Output Parameters: 414c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4153c632250SBarry Smith . w - new iterate 41629e0b56fSBarry Smith . gnorm - not changed 41729e0b56fSBarry Smith . ynorm - not changed 418a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 419fee21e36SBarry Smith 42029e0b56fSBarry Smith Options Database Key: 42117bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 42229e0b56fSBarry Smith 4238cbba510SBarry Smith Notes: 42417bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 425ea56f5baSLois Curfman McInnes either the options 426ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 427ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 428329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 429329f5518SBarry Smith otherwise, the SNES solver will generate an error. 430329f5518SBarry Smith 431329f5518SBarry Smith During the final iteration this will not evaluate the function at 432329f5518SBarry Smith the solution point. This is to save a function evaluation while 433329f5518SBarry Smith using pseudo-timestepping. 4348cbba510SBarry Smith 435ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 436a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 437ea56f5baSLois Curfman McInnes correct, since they are not computed. 438ea56f5baSLois Curfman McInnes 439ea56f5baSLois Curfman McInnes Level: advanced 4408cbba510SBarry Smith 44129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 44229e0b56fSBarry Smith 44317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 44417bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 44529e0b56fSBarry Smith @*/ 446dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 44729e0b56fSBarry Smith { 448dfbe8321SBarry Smith PetscErrorCode ierr; 4494b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4503c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 45129e0b56fSBarry Smith 4523a40ed3dSBarry Smith PetscFunctionBegin; 453a7cc72afSBarry Smith *flag = PETSC_TRUE; 454d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 45579f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4563c632250SBarry Smith if (neP->postcheckstep) { 4573c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4583c632250SBarry Smith } 4593c632250SBarry Smith if (changed_y) { 46079f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4618f99978dSLois Curfman McInnes } 462329f5518SBarry Smith 463329f5518SBarry Smith /* don't evaluate function the last time through */ 464329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4654936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 466329f5518SBarry Smith } 467d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 46929e0b56fSBarry Smith } 47004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4714a2ae208SSatish Balay #undef __FUNCT__ 47217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4734b828684SBarry Smith /*@C 47417bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4755e76c431SBarry Smith 476c7afd0dbSLois Curfman McInnes Collective on SNES 477c7afd0dbSLois Curfman McInnes 4785e76c431SBarry Smith Input Parameters: 479c7afd0dbSLois Curfman McInnes + snes - nonlinear context 480af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4815e76c431SBarry Smith . x - current iterate 4825e76c431SBarry Smith . f - residual evaluated at x 4833c632250SBarry Smith . y - search direction 4845e76c431SBarry Smith . w - work vector 485dc357ad2SBarry Smith . fnorm - 2-norm of f 486dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 4875e76c431SBarry Smith 488393d2d9aSLois Curfman McInnes Output Parameters: 489c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4903c632250SBarry Smith . w - new iterate 4915e76c431SBarry Smith . gnorm - 2-norm of g 4925e76c431SBarry Smith . ynorm - 2-norm of search length 493a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 494fee21e36SBarry Smith 495c4a48953SLois Curfman McInnes Options Database Key: 49617bae607SBarry Smith $ -snes_ls cubic - Activates SNESLineSearchCubic() 497c4a48953SLois Curfman McInnes 4985e76c431SBarry Smith Notes: 4995e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 5005e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 5015e76c431SBarry Smith 50236851e7fSLois Curfman McInnes Level: advanced 50336851e7fSLois Curfman McInnes 50428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 50528ae5a21SLois Curfman McInnes 50617bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 5075e76c431SBarry Smith @*/ 508dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 5095e76c431SBarry Smith { 510406556e6SLois Curfman McInnes /* 511406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 512406556e6SLois Curfman McInnes minimization problem: 513406556e6SLois Curfman McInnes min z(x): R^n -> R, 514406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 515406556e6SLois Curfman McInnes */ 516406556e6SLois Curfman McInnes 5175ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 518275d25c3SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp; 519aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 520f5b6597dSBarry Smith PetscScalar cinitslope; 5216b5873e3SBarry Smith #endif 522dfbe8321SBarry Smith PetscErrorCode ierr; 523a7cc72afSBarry Smith PetscInt count; 5244b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 52579f36460SBarry Smith PetscScalar scale; 5263c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 5275e76c431SBarry Smith 5283a40ed3dSBarry Smith PetscFunctionBegin; 529d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 530a7cc72afSBarry Smith *flag = PETSC_TRUE; 5315e76c431SBarry Smith alpha = neP->alpha; 5325e76c431SBarry Smith maxstep = neP->maxstep; 533dc357ad2SBarry Smith steptol = snes->xtol; 5345e76c431SBarry Smith 535cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 5369e247f21SBarry Smith if (!*ynorm) { 537ae15b995SBarry Smith ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 538a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5393c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 540a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 541a7cc72afSBarry Smith *flag = PETSC_FALSE; 542ad922ac9SBarry Smith goto theend1; 54394a9d846SBarry Smith } 5445e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5455e42470aSBarry Smith scale = maxstep/(*ynorm); 546ae15b995SBarry Smith ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr); 5472dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 5485e76c431SBarry Smith *ynorm = maxstep; 5495e76c431SBarry Smith } 550ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5515ca10a37SBarry Smith minlambda = steptol/rellength; 552a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 553aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 554a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 555329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5565e42470aSBarry Smith #else 557a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5585e42470aSBarry Smith #endif 5595e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5605e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5615e76c431SBarry Smith 562e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 56343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 564ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 56543e71028SBarry Smith *flag = PETSC_FALSE; 56643e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 56743e71028SBarry Smith goto theend1; 56843e71028SBarry Smith } 5694936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5704936397dSBarry Smith if (snes->domainerror) { 5714936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5724936397dSBarry Smith PetscFunctionReturn(0); 57319717074SBarry Smith } 574313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 575f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 576f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 5775d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 578f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 579ad922ac9SBarry Smith goto theend1; 5805e76c431SBarry Smith } 581dc357ad2SBarry Smith if (*ynorm < xnorm*snes->xtol) { 582dc357ad2SBarry Smith ierr = PetscInfo3(snes,"Using full step: because ynorm %G < xnorm %G * steptol %G (i.e. Newton step is below tolerance)\n",*ynorm,xnorm,snes->xtol);CHKERRQ(ierr); 583dc357ad2SBarry Smith *flag = PETSC_TRUE; 584dc357ad2SBarry Smith goto theend1; 585dc357ad2SBarry Smith } 586dc357ad2SBarry Smith 5875e76c431SBarry Smith 5885e76c431SBarry Smith /* Fit points with quadratic */ 589313b5538SBarry Smith lambda = 1.0; 590a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5915e76c431SBarry Smith lambdaprev = lambda; 5925e76c431SBarry Smith gnormprev = *gnorm; 593329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 594ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 595ddd12767SBarry Smith else lambda = lambdatemp; 596275d25c3SBarry Smith 597e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 59843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 599ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 60043e71028SBarry Smith *flag = PETSC_FALSE; 60143e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 60243e71028SBarry Smith goto theend1; 60343e71028SBarry Smith } 604f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6054936397dSBarry Smith if (snes->domainerror) { 6064936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6074936397dSBarry Smith PetscFunctionReturn(0); 60819717074SBarry Smith } 609cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 610f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 611f5b6597dSBarry Smith ierr = PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 6125ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 613ae15b995SBarry Smith ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 614ad922ac9SBarry Smith goto theend1; 6155e76c431SBarry Smith } 6165e76c431SBarry Smith 6175e76c431SBarry Smith /* Fit points with cubic */ 6185e76c431SBarry Smith count = 1; 619*bb9195b6SBarry Smith while (PETSC_TRUE) { 620dc357ad2SBarry Smith if (lambda <= minlambda) { 621dc357ad2SBarry Smith ierr = PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 622ae15b995SBarry Smith ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 62343e71028SBarry Smith *flag = PETSC_FALSE; 62443e71028SBarry Smith break; 6255e76c431SBarry Smith } 6265d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6275d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 628ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6292b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6305e76c431SBarry Smith d = b*b - 3*a*initslope; 6315e76c431SBarry Smith if (d < 0.0) d = 0.0; 6325e76c431SBarry Smith if (a == 0.0) { 6335e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6345e76c431SBarry Smith } else { 6355e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6365e76c431SBarry Smith } 6375e76c431SBarry Smith lambdaprev = lambda; 6385e76c431SBarry Smith gnormprev = *gnorm; 639329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 640329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6415e76c431SBarry Smith else lambda = lambdatemp; 642e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 64343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 644ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 645ae15b995SBarry Smith ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 646ed98166eSMatthew Knepley *flag = PETSC_FALSE; 64743e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 648ed98166eSMatthew Knepley break; 649ed98166eSMatthew Knepley } 650f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6514936397dSBarry Smith if (snes->domainerror) { 6524936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6534936397dSBarry Smith PetscFunctionReturn(0); 65419717074SBarry Smith } 655cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 656f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6575ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 658f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 65943e71028SBarry Smith break; 6602b022350SLois Curfman McInnes } else { 661f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 6625e76c431SBarry Smith } 6635e76c431SBarry Smith count++; 6645e76c431SBarry Smith } 665ad922ac9SBarry Smith theend1: 6668f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6673c632250SBarry Smith if (neP->postcheckstep && *flag) { 6683c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6693c632250SBarry Smith if (changed_y) { 67079f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6713c632250SBarry Smith } 6723c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 673f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6744936397dSBarry Smith if (snes->domainerror) { 6754936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6764936397dSBarry Smith PetscFunctionReturn(0); 67719717074SBarry Smith } 6788f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 679f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 680a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6818f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 682a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6838f99978dSLois Curfman McInnes } 6848f99978dSLois Curfman McInnes } 685d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6863a40ed3dSBarry Smith PetscFunctionReturn(0); 6875e76c431SBarry Smith } 68804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6894a2ae208SSatish Balay #undef __FUNCT__ 69017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 6914b828684SBarry Smith /*@C 69217bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 6935e76c431SBarry Smith 694c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 695c7afd0dbSLois Curfman McInnes 6965e42470aSBarry Smith Input Parameters: 697c7afd0dbSLois Curfman McInnes + snes - the SNES context 698af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6995e42470aSBarry Smith . x - current iterate 7005e42470aSBarry Smith . f - residual evaluated at x 7013c632250SBarry Smith . y - search direction 7025e42470aSBarry Smith . w - work vector 703dc357ad2SBarry Smith . fnorm - 2-norm of f 704dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 7055e42470aSBarry Smith 706c4a48953SLois Curfman McInnes Output Parameters: 7077f3332b4SBarry Smith + g - residual evaluated at new iterate w 7087f3332b4SBarry Smith . w - new iterate (x + alpha*y) 7095e42470aSBarry Smith . gnorm - 2-norm of g 7105e42470aSBarry Smith . ynorm - 2-norm of search length 711a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 712fee21e36SBarry Smith 713ce9499c7SBarry Smith Options Database Keys: 714ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 715ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 716ce9499c7SBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 717dc357ad2SBarry Smith - -snes_stol <steptol> - Sets steptol, this is the minimum step size that the line search code 718ce9499c7SBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 719ce9499c7SBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 7205e42470aSBarry Smith Notes: 72117bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7225e42470aSBarry Smith 72336851e7fSLois Curfman McInnes Level: advanced 72436851e7fSLois Curfman McInnes 725f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 726f59ffdedSLois Curfman McInnes 72717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7285e42470aSBarry Smith @*/ 729dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 7305e76c431SBarry Smith { 731406556e6SLois Curfman McInnes /* 732406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 733406556e6SLois Curfman McInnes minimization problem: 734406556e6SLois Curfman McInnes min z(x): R^n -> R, 735406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 736406556e6SLois Curfman McInnes */ 737275d25c3SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength; 738aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 739f5b6597dSBarry Smith PetscScalar cinitslope; 7406b5873e3SBarry Smith #endif 741dfbe8321SBarry Smith PetscErrorCode ierr; 742a7cc72afSBarry Smith PetscInt count; 7434b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 74479f36460SBarry Smith PetscScalar scale; 7453c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7465e76c431SBarry Smith 7473a40ed3dSBarry Smith PetscFunctionBegin; 748d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 749a7cc72afSBarry Smith *flag = PETSC_TRUE; 7505e42470aSBarry Smith alpha = neP->alpha; 7515e42470aSBarry Smith maxstep = neP->maxstep; 752dc357ad2SBarry Smith steptol = snes->xtol; 7535e76c431SBarry Smith 7543505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 75563b9fa5eSBarry Smith if (*ynorm == 0.0) { 756ae15b995SBarry Smith ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 757b37302c6SLois Curfman McInnes *gnorm = fnorm; 758e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 759b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 760a7cc72afSBarry Smith *flag = PETSC_FALSE; 761ad922ac9SBarry Smith goto theend2; 76294a9d846SBarry Smith } 7635e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 7645e42470aSBarry Smith scale = maxstep/(*ynorm); 7652dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 7665e42470aSBarry Smith *ynorm = maxstep; 7675e76c431SBarry Smith } 768ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 7695ca10a37SBarry Smith minlambda = steptol/rellength; 770a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 771aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 772a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 773329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7745e42470aSBarry Smith #else 775a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7765e42470aSBarry Smith #endif 7775e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7785e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 779f8833479SBarry Smith ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 7805e42470aSBarry Smith 781e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 78243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 783ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 78443e71028SBarry Smith *flag = PETSC_FALSE; 78543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 78643e71028SBarry Smith goto theend2; 78743e71028SBarry Smith } 7884936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 7894936397dSBarry Smith if (snes->domainerror) { 7904936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7914936397dSBarry Smith PetscFunctionReturn(0); 79219717074SBarry Smith } 793cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 794f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 7955d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 796ae15b995SBarry Smith ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 797ad922ac9SBarry Smith goto theend2; 7985e42470aSBarry Smith } 799dc357ad2SBarry Smith if (*ynorm < xnorm*snes->xtol) { 800dc357ad2SBarry Smith ierr = PetscInfo3(snes,"Using full step: because ynorm %G < xnorm %G * steptol %G (i.e. Newton step is below tolerance)\n",*ynorm,xnorm,snes->xtol);CHKERRQ(ierr); 801dc357ad2SBarry Smith *flag = PETSC_TRUE; 802dc357ad2SBarry Smith goto theend2; 803dc357ad2SBarry Smith } 8045e42470aSBarry Smith 8055e42470aSBarry Smith /* Fit points with quadratic */ 806313b5538SBarry Smith lambda = 1.0; 8075e42470aSBarry Smith count = 1; 8085ca10a37SBarry Smith while (PETSC_TRUE) { 8095e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 810ae15b995SBarry Smith ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 811ae15b995SBarry Smith ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 812e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 81343e71028SBarry Smith *flag = PETSC_FALSE; 81443e71028SBarry Smith break; 8155e42470aSBarry Smith } 816a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 817329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 818329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 819329f5518SBarry Smith else lambda = lambdatemp; 820275d25c3SBarry Smith 821e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 82243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 823ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 824ae15b995SBarry Smith ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 825ed98166eSMatthew Knepley *flag = PETSC_FALSE; 82643e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 827ed98166eSMatthew Knepley break; 828ed98166eSMatthew Knepley } 829f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8304936397dSBarry Smith if (snes->domainerror) { 8314936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8324936397dSBarry Smith PetscFunctionReturn(0); 83319717074SBarry Smith } 834cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 835f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8365ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 837ae15b995SBarry Smith ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 83806259719SBarry Smith break; 8395e42470aSBarry Smith } 8405e42470aSBarry Smith count++; 8415e42470aSBarry Smith } 842ad922ac9SBarry Smith theend2: 8438f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8443c632250SBarry Smith if (neP->postcheckstep) { 8453c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8463c632250SBarry Smith if (changed_y) { 84779f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8483c632250SBarry Smith } 8493c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8503c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 8514936397dSBarry Smith if (snes->domainerror) { 8524936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8534936397dSBarry Smith PetscFunctionReturn(0); 85419717074SBarry Smith } 8558f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 856a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8578f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 858a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 859f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8608f99978dSLois Curfman McInnes } 8618f99978dSLois Curfman McInnes } 862d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8633a40ed3dSBarry Smith PetscFunctionReturn(0); 8645e76c431SBarry Smith } 8652343ba6eSBarry Smith 86604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8674a2ae208SSatish Balay #undef __FUNCT__ 86817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 869c9e6a524SLois Curfman McInnes /*@C 87017bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8714b27c08aSLois Curfman McInnes by the method SNESLS. 8725e76c431SBarry Smith 8735e76c431SBarry Smith Input Parameters: 874c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 875af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 876c7afd0dbSLois Curfman McInnes - func - pointer to int function 8775e76c431SBarry Smith 878fee21e36SBarry Smith Collective on SNES 879fee21e36SBarry Smith 880c4a48953SLois Curfman McInnes Available Routines: 88117bae607SBarry Smith + SNESLineSearchCubic() - default line search 88217bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 88317bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 88417bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8855e76c431SBarry Smith 886c4a48953SLois Curfman McInnes Options Database Keys: 8874b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8884b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 8894b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 890dc357ad2SBarry Smith - -snes_xtol <steptol> - Sets steptol, this is the minimum step size that the line search code. This is the same as the minimum step length 8913304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 892c4a48953SLois Curfman McInnes 8935e76c431SBarry Smith Calling sequence of func: 894bd208895SLois Curfman McInnes .vb 895dc357ad2SBarry Smith func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 896bd208895SLois Curfman McInnes .ve 8975e76c431SBarry Smith 8985e76c431SBarry Smith Input parameters for func: 899c7afd0dbSLois Curfman McInnes + snes - nonlinear context 900af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 9015e76c431SBarry Smith . x - current iterate 9025e76c431SBarry Smith . f - residual evaluated at x 9033c632250SBarry Smith . y - search direction 9045e76c431SBarry Smith . w - work vector 905c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9065e76c431SBarry Smith 9075e76c431SBarry Smith Output parameters for func: 908c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9093c632250SBarry Smith . w - new iterate 9105e76c431SBarry Smith . gnorm - 2-norm of g 9115e76c431SBarry Smith . ynorm - 2-norm of search length 912a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 913f59ffdedSLois Curfman McInnes 91436851e7fSLois Curfman McInnes Level: advanced 91536851e7fSLois Curfman McInnes 916f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 917f59ffdedSLois Curfman McInnes 91817bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 91917bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9205e76c431SBarry Smith @*/ 921dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 9225e76c431SBarry Smith { 923dc357ad2SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 92482bf6240SBarry Smith 9253a40ed3dSBarry Smith PetscFunctionBegin; 92617bae607SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 92782bf6240SBarry Smith if (f) { 928af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 92982bf6240SBarry Smith } 9303a40ed3dSBarry Smith PetscFunctionReturn(0); 9315e76c431SBarry Smith } 9328e019c35SBarry Smith 933dc357ad2SBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 93404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 935fb2e594dSBarry Smith EXTERN_C_BEGIN 9364a2ae208SSatish Balay #undef __FUNCT__ 93717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 93817bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 93982bf6240SBarry Smith { 94082bf6240SBarry Smith PetscFunctionBegin; 9414b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9424b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 94382bf6240SBarry Smith PetscFunctionReturn(0); 94482bf6240SBarry Smith } 945fb2e594dSBarry Smith EXTERN_C_END 94604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9474a2ae208SSatish Balay #undef __FUNCT__ 9483c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 949c8dd0c0dSLois Curfman McInnes /*@C 9503c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9514b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 952c8dd0c0dSLois Curfman McInnes 953c8dd0c0dSLois Curfman McInnes Input Parameters: 954c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9553c632250SBarry Smith . func - pointer to function 956c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 957c8dd0c0dSLois Curfman McInnes 958c8dd0c0dSLois Curfman McInnes Collective on SNES 959c8dd0c0dSLois Curfman McInnes 960c8dd0c0dSLois Curfman McInnes Calling sequence of func: 961c8dd0c0dSLois Curfman McInnes .vb 9623c632250SBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 963c8dd0c0dSLois Curfman McInnes .ve 964b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 965b1ae27eaSLois Curfman McInnes on failure. 966c8dd0c0dSLois Curfman McInnes 967c8dd0c0dSLois Curfman McInnes Input parameters for func: 968c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 969c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9703c632250SBarry Smith . x - previous iterate 9713c632250SBarry Smith . y - new search direction and length 9723c632250SBarry Smith - w - current candidate iterate 973c8dd0c0dSLois Curfman McInnes 974c8dd0c0dSLois Curfman McInnes Output parameters for func: 9753c632250SBarry Smith + y - search direction (possibly changed) 9763c632250SBarry Smith . w - current iterate (possibly modified) 9773c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 9783c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 979c8dd0c0dSLois Curfman McInnes 980c8dd0c0dSLois Curfman McInnes Level: advanced 981c8dd0c0dSLois Curfman McInnes 9829e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 9839e247f21SBarry Smith 9843c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 9853c632250SBarry Smith 9863c632250SBarry Smith On input w = x + y 9873c632250SBarry Smith 98817bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 989b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 990ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9918f99978dSLois Curfman McInnes 99217bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 993ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 994ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 995ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 9969e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 997b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9988f99978dSLois Curfman McInnes 999c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 1000c8dd0c0dSLois Curfman McInnes 100117bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1002c8dd0c0dSLois Curfman McInnes @*/ 10033c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 1004c8dd0c0dSLois Curfman McInnes { 10053c632250SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 1006c8dd0c0dSLois Curfman McInnes 1007c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10083c632250SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1009c8dd0c0dSLois Curfman McInnes if (f) { 1010c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1011c8dd0c0dSLois Curfman McInnes } 1012c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1013c8dd0c0dSLois Curfman McInnes } 10149c3ca13aSBarry Smith #undef __FUNCT__ 10159c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10169c3ca13aSBarry Smith /*@C 10179c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10187e4bb74cSBarry Smith before the line search is called. 10199c3ca13aSBarry Smith 10209c3ca13aSBarry Smith Input Parameters: 10219c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10229c3ca13aSBarry Smith . func - pointer to function 10239c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10249c3ca13aSBarry Smith 10259c3ca13aSBarry Smith Collective on SNES 10269c3ca13aSBarry Smith 10279c3ca13aSBarry Smith Calling sequence of func: 10289c3ca13aSBarry Smith .vb 10291e3347e8SBarry Smith int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 10309c3ca13aSBarry Smith .ve 10319c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10329c3ca13aSBarry Smith on failure. 10339c3ca13aSBarry Smith 10349c3ca13aSBarry Smith Input parameters for func: 10359c3ca13aSBarry Smith + snes - nonlinear context 10369c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10379c3ca13aSBarry Smith . x - previous iterate 10389c3ca13aSBarry Smith - y - new search direction and length 10399c3ca13aSBarry Smith 10409c3ca13aSBarry Smith Output parameters for func: 10419c3ca13aSBarry Smith + y - search direction (possibly changed) 10429c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10439c3ca13aSBarry Smith 10449c3ca13aSBarry Smith Level: advanced 10459c3ca13aSBarry Smith 10469c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10479c3ca13aSBarry Smith 10489c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10499c3ca13aSBarry Smith 10507e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 10519c3ca13aSBarry Smith @*/ 10529c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 10539c3ca13aSBarry Smith { 10549c3ca13aSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 10559c3ca13aSBarry Smith 10569c3ca13aSBarry Smith PetscFunctionBegin; 10579c3ca13aSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 10589c3ca13aSBarry Smith if (f) { 10599c3ca13aSBarry Smith ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 10609c3ca13aSBarry Smith } 10619c3ca13aSBarry Smith PetscFunctionReturn(0); 10629c3ca13aSBarry Smith } 10639c3ca13aSBarry Smith 1064c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 10659c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 10669c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1067c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 10684a2ae208SSatish Balay #undef __FUNCT__ 10693c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 10709c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1071c8dd0c0dSLois Curfman McInnes { 1072c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10733c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 10743c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1075c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1076c8dd0c0dSLois Curfman McInnes } 1077c8dd0c0dSLois Curfman McInnes EXTERN_C_END 10789c3ca13aSBarry Smith 10799c3ca13aSBarry Smith EXTERN_C_BEGIN 10809c3ca13aSBarry Smith #undef __FUNCT__ 10819c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 10829c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 10839c3ca13aSBarry Smith { 10849c3ca13aSBarry Smith PetscFunctionBegin; 10859c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 10869c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 10879c3ca13aSBarry Smith PetscFunctionReturn(0); 10889c3ca13aSBarry Smith } 10899c3ca13aSBarry Smith EXTERN_C_END 1090329e7e40SMatthew Knepley 1091329e7e40SMatthew Knepley /* 10924b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 109304d965bbSLois Curfman McInnes 109404d965bbSLois Curfman McInnes Input Parameters: 109504d965bbSLois Curfman McInnes . SNES - the SNES context 109604d965bbSLois Curfman McInnes . viewer - visualization context 109704d965bbSLois Curfman McInnes 109804d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 109904d965bbSLois Curfman McInnes */ 11004a2ae208SSatish Balay #undef __FUNCT__ 11014b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11026849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1103a935fc98SLois Curfman McInnes { 11044b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11052fc52814SBarry Smith const char *cstr; 1106dfbe8321SBarry Smith PetscErrorCode ierr; 110732077d6dSBarry Smith PetscTruth iascii; 1108a935fc98SLois Curfman McInnes 11093a40ed3dSBarry Smith PetscFunctionBegin; 111032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 111132077d6dSBarry Smith if (iascii) { 111217bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 111317bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 111417bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 111519bcc07fSBarry Smith else cstr = "unknown"; 1116b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1117dc357ad2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G\n",ls->alpha,ls->maxstep);CHKERRQ(ierr); 11185cd90555SBarry Smith } else { 111979a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 112019bcc07fSBarry Smith } 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 1122a935fc98SLois Curfman McInnes } 112304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 112404d965bbSLois Curfman McInnes /* 11254b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 112604d965bbSLois Curfman McInnes 112704d965bbSLois Curfman McInnes Input Parameter: 112804d965bbSLois Curfman McInnes . snes - the SNES context 112904d965bbSLois Curfman McInnes 113004d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 113104d965bbSLois Curfman McInnes */ 11324a2ae208SSatish Balay #undef __FUNCT__ 11334b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11346849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11355e42470aSBarry Smith { 11364b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1137e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1138dfbe8321SBarry Smith PetscErrorCode ierr; 1139a7cc72afSBarry Smith PetscInt indx; 1140f1af5d2fSBarry Smith PetscTruth flg; 11415e42470aSBarry Smith 11423a40ed3dSBarry Smith PetscFunctionBegin; 1143b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 11444b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 11454b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1146186905e3SBarry Smith 114717bae607SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 114825cdf11fSBarry Smith if (flg) { 114922e36657SBarry Smith switch (indx) { 1150b49fd9e1SBarry Smith case 0: 115117bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1152b49fd9e1SBarry Smith break; 1153b49fd9e1SBarry Smith case 1: 115417bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1155b49fd9e1SBarry Smith break; 1156b49fd9e1SBarry Smith case 2: 115717bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1158b49fd9e1SBarry Smith break; 1159b49fd9e1SBarry Smith case 3: 116017bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1161b49fd9e1SBarry Smith break; 11625e42470aSBarry Smith } 11635e42470aSBarry Smith } 1164b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11653a40ed3dSBarry Smith PetscFunctionReturn(0); 11665e42470aSBarry Smith } 116704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1168ebe3b25bSBarry Smith /*MC 1169ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 117004d965bbSLois Curfman McInnes 1171ebe3b25bSBarry Smith Options Database: 1172ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1173ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1174ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1175ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1176ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1177ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 117804d965bbSLois Curfman McInnes 1179ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 118004d965bbSLois Curfman McInnes 1181ee3001cbSBarry Smith Level: beginner 1182ee3001cbSBarry Smith 118317bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 118417bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 118517bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1186ebe3b25bSBarry Smith 1187ebe3b25bSBarry Smith M*/ 1188fb2e594dSBarry Smith EXTERN_C_BEGIN 11894a2ae208SSatish Balay #undef __FUNCT__ 11904b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 119163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 11925e42470aSBarry Smith { 1193dfbe8321SBarry Smith PetscErrorCode ierr; 11944b27c08aSLois Curfman McInnes SNES_LS *neP; 11955e42470aSBarry Smith 11963a40ed3dSBarry Smith PetscFunctionBegin; 1197e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 1198e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 1199e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 1200e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 1201e7788613SBarry Smith snes->ops->view = SNESView_LS; 12025e42470aSBarry Smith 120338f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 12045e42470aSBarry Smith snes->data = (void*)neP; 12055e42470aSBarry Smith neP->alpha = 1.e-4; 12065e42470aSBarry Smith neP->maxstep = 1.e8; 120717bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1208c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12093c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12103c632250SBarry Smith neP->postcheck = PETSC_NULL; 12113c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12123c632250SBarry Smith neP->precheck = PETSC_NULL; 121382bf6240SBarry Smith 1214557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1215557d3f75SLisandro Dalcin "SNESLineSearchSet_LS", 1216557d3f75SLisandro Dalcin SNESLineSearchSet_LS);CHKERRQ(ierr); 1217557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1218557d3f75SLisandro Dalcin "SNESLineSearchSetPostCheck_LS", 1219557d3f75SLisandro Dalcin SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1220557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1221557d3f75SLisandro Dalcin "SNESLineSearchSetPreCheck_LS", 1222557d3f75SLisandro Dalcin SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 122382bf6240SBarry Smith 12243a40ed3dSBarry Smith PetscFunctionReturn(0); 12255e42470aSBarry Smith } 1226fb2e594dSBarry Smith EXTERN_C_END 122782bf6240SBarry Smith 122882bf6240SBarry Smith 122982bf6240SBarry Smith 1230