163dd3a1aSKris Buschelman #define PETSCSNES_DLL 25e76c431SBarry Smith 3724775ccSBarry Smith #include "../src/snes/impls/ls/lsimpl.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); 6775567043SBarry Smith if (a1 != 0.0) { 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; 15775567043SBarry Smith snes->norm = 0.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 { 30594298653SBarry Smith SNES_LS *ls = (SNES_LS*) snes->data; 306dfbe8321SBarry Smith PetscErrorCode ierr; 3073a40ed3dSBarry Smith 3083a40ed3dSBarry Smith PetscFunctionBegin; 30985385478SLisandro Dalcin if (snes->vec_sol_update) { 31085385478SLisandro Dalcin ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 31185385478SLisandro Dalcin snes->vec_sol_update = PETSC_NULL; 31285385478SLisandro Dalcin } 3135baf8537SBarry Smith if (snes->nwork) { 3144b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 31585385478SLisandro Dalcin snes->nwork = 0; 31685385478SLisandro Dalcin snes->work = PETSC_NULL; 31721c89e3eSBarry Smith } 31894298653SBarry Smith if (ls->monitor) { 31994298653SBarry Smith ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 32094298653SBarry Smith } 3215c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 322557d3f75SLisandro Dalcin 323557d3f75SLisandro Dalcin /* clear composed functions */ 324557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 325557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 326557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 327557d3f75SLisandro Dalcin 3283a40ed3dSBarry Smith PetscFunctionReturn(0); 3295e76c431SBarry Smith } 33004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3314a2ae208SSatish Balay #undef __FUNCT__ 33217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 33382bf6240SBarry Smith 3344b828684SBarry Smith /*@C 33517bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3365e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3375e76c431SBarry Smith to serve as a template and is not recommended for general use. 3385e76c431SBarry Smith 339c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 340c7afd0dbSLois Curfman McInnes 3415e76c431SBarry Smith Input Parameters: 342c7afd0dbSLois Curfman McInnes + snes - nonlinear context 343af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3445e76c431SBarry Smith . x - current iterate 3455e76c431SBarry Smith . f - residual evaluated at x 3463c632250SBarry Smith . y - search direction 347dc357ad2SBarry Smith . fnorm - 2-norm of f 348dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 3495e76c431SBarry Smith 350c4a48953SLois Curfman McInnes Output Parameters: 351c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3523c632250SBarry Smith . w - new iterate 3535e76c431SBarry Smith . gnorm - 2-norm of g 3545e76c431SBarry Smith . ynorm - 2-norm of search length 355a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 356fee21e36SBarry Smith 357c4a48953SLois Curfman McInnes Options Database Key: 35817bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 359c4a48953SLois Curfman McInnes 36036851e7fSLois Curfman McInnes Level: advanced 36136851e7fSLois Curfman McInnes 36228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 36328ae5a21SLois Curfman McInnes 36417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 36517bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3665e76c431SBarry Smith @*/ 367dc357ad2SBarry 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) 3685e76c431SBarry Smith { 369dfbe8321SBarry Smith PetscErrorCode ierr; 3704b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3713c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3725334005bSBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 374a7cc72afSBarry Smith *flag = PETSC_TRUE; 375d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 376a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 37779f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3783c632250SBarry Smith if (neP->postcheckstep) { 3793c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3808f99978dSLois Curfman McInnes } 3813c632250SBarry Smith if (changed_y) { 38279f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3833c632250SBarry Smith } 3844936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 3854936397dSBarry Smith if (!snes->domainerror) { 386a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 387f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3884936397dSBarry Smith } 389d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3903a40ed3dSBarry Smith PetscFunctionReturn(0); 3915e76c431SBarry Smith } 39204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3934a2ae208SSatish Balay #undef __FUNCT__ 39417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 39582bf6240SBarry Smith 39629e0b56fSBarry Smith /*@C 39717bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 39829e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 39929e0b56fSBarry Smith even compute the norm of the function or search direction; this 40029e0b56fSBarry Smith is intended only when you know the full step is fine and are 40129e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 402c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 403c7afd0dbSLois Curfman McInnes 404c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 40529e0b56fSBarry Smith 40629e0b56fSBarry Smith Input Parameters: 407c7afd0dbSLois Curfman McInnes + snes - nonlinear context 408af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 40929e0b56fSBarry Smith . x - current iterate 41029e0b56fSBarry Smith . f - residual evaluated at x 4113c632250SBarry Smith . y - search direction 41229e0b56fSBarry Smith . w - work vector 413dc357ad2SBarry Smith . fnorm - 2-norm of f 414dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 41529e0b56fSBarry Smith 41629e0b56fSBarry Smith Output Parameters: 417c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4183c632250SBarry Smith . w - new iterate 41929e0b56fSBarry Smith . gnorm - not changed 42029e0b56fSBarry Smith . ynorm - not changed 421a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 422fee21e36SBarry Smith 42329e0b56fSBarry Smith Options Database Key: 42417bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 42529e0b56fSBarry Smith 4268cbba510SBarry Smith Notes: 42717bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 428ea56f5baSLois Curfman McInnes either the options 429ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 430ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 431329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 432329f5518SBarry Smith otherwise, the SNES solver will generate an error. 433329f5518SBarry Smith 434329f5518SBarry Smith During the final iteration this will not evaluate the function at 435329f5518SBarry Smith the solution point. This is to save a function evaluation while 436329f5518SBarry Smith using pseudo-timestepping. 4378cbba510SBarry Smith 438ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 439a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 440ea56f5baSLois Curfman McInnes correct, since they are not computed. 441ea56f5baSLois Curfman McInnes 442ea56f5baSLois Curfman McInnes Level: advanced 4438cbba510SBarry Smith 44429e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 44529e0b56fSBarry Smith 44617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 44717bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 44829e0b56fSBarry Smith @*/ 449dc357ad2SBarry 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) 45029e0b56fSBarry Smith { 451dfbe8321SBarry Smith PetscErrorCode ierr; 4524b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4533c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 45429e0b56fSBarry Smith 4553a40ed3dSBarry Smith PetscFunctionBegin; 456a7cc72afSBarry Smith *flag = PETSC_TRUE; 457d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 45879f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4593c632250SBarry Smith if (neP->postcheckstep) { 4603c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4613c632250SBarry Smith } 4623c632250SBarry Smith if (changed_y) { 46379f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4648f99978dSLois Curfman McInnes } 465329f5518SBarry Smith 466329f5518SBarry Smith /* don't evaluate function the last time through */ 467329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4684936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 469329f5518SBarry Smith } 470d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4713a40ed3dSBarry Smith PetscFunctionReturn(0); 47229e0b56fSBarry Smith } 47304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4744a2ae208SSatish Balay #undef __FUNCT__ 47517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4764b828684SBarry Smith /*@C 47717bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4785e76c431SBarry Smith 479c7afd0dbSLois Curfman McInnes Collective on SNES 480c7afd0dbSLois Curfman McInnes 4815e76c431SBarry Smith Input Parameters: 482c7afd0dbSLois Curfman McInnes + snes - nonlinear context 483af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4845e76c431SBarry Smith . x - current iterate 4855e76c431SBarry Smith . f - residual evaluated at x 4863c632250SBarry Smith . y - search direction 4875e76c431SBarry Smith . w - work vector 488dc357ad2SBarry Smith . fnorm - 2-norm of f 489dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 4905e76c431SBarry Smith 491393d2d9aSLois Curfman McInnes Output Parameters: 492c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4933c632250SBarry Smith . w - new iterate 4945e76c431SBarry Smith . gnorm - 2-norm of g 4955e76c431SBarry Smith . ynorm - 2-norm of search length 496a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 497fee21e36SBarry Smith 498c4a48953SLois Curfman McInnes Options Database Key: 499e106eecfSBarry Smith + -snes_ls cubic - Activates SNESLineSearchCubic() 500e106eecfSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 501e106eecfSBarry Smith . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 502e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 503e106eecfSBarry Smith 504c4a48953SLois Curfman McInnes 5055e76c431SBarry Smith Notes: 5065e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 5075e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 5085e76c431SBarry Smith 50936851e7fSLois Curfman McInnes Level: advanced 51036851e7fSLois Curfman McInnes 51128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 51228ae5a21SLois Curfman McInnes 51317bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 5145e76c431SBarry Smith @*/ 515dc357ad2SBarry 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) 5165e76c431SBarry Smith { 517406556e6SLois Curfman McInnes /* 518406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 519406556e6SLois Curfman McInnes minimization problem: 520406556e6SLois Curfman McInnes min z(x): R^n -> R, 521406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 522406556e6SLois Curfman McInnes */ 523406556e6SLois Curfman McInnes 524e106eecfSBarry Smith PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 525e106eecfSBarry Smith PetscReal minlambda,lambda,lambdatemp; 526aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 527f5b6597dSBarry Smith PetscScalar cinitslope; 5286b5873e3SBarry Smith #endif 529dfbe8321SBarry Smith PetscErrorCode ierr; 530a7cc72afSBarry Smith PetscInt count; 5314b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 5323c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 53394298653SBarry Smith MPI_Comm comm; 5345e76c431SBarry Smith 5353a40ed3dSBarry Smith PetscFunctionBegin; 53694298653SBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 537d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 538a7cc72afSBarry Smith *flag = PETSC_TRUE; 5395e76c431SBarry Smith 540cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 54175567043SBarry Smith if (*ynorm == 0.0) { 54294298653SBarry Smith if (neP->monitor) { 54394298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 54494298653SBarry Smith } 545a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5463c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 547a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 548a7cc72afSBarry Smith *flag = PETSC_FALSE; 549ad922ac9SBarry Smith goto theend1; 55094a9d846SBarry Smith } 551e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 55294298653SBarry Smith if (neP->monitor) { 55394298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 55494298653SBarry Smith } 555e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 556e106eecfSBarry Smith *ynorm = neP->maxstep; 5575e76c431SBarry Smith } 558ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 559e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 560a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 561aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 562a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 563329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5645e42470aSBarry Smith #else 565a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5665e42470aSBarry Smith #endif 5675e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5685e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5695e76c431SBarry Smith 570e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 57143e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 572ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 57343e71028SBarry Smith *flag = PETSC_FALSE; 57443e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 57543e71028SBarry Smith goto theend1; 57643e71028SBarry Smith } 5774936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5784936397dSBarry Smith if (snes->domainerror) { 5794936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5804936397dSBarry Smith PetscFunctionReturn(0); 58119717074SBarry Smith } 582313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 583f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 584f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 585e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 58694298653SBarry Smith if (neP->monitor) { 58794298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 58894298653SBarry Smith } 589ad922ac9SBarry Smith goto theend1; 5905e76c431SBarry Smith } 5915e76c431SBarry Smith 5925e76c431SBarry Smith /* Fit points with quadratic */ 593313b5538SBarry Smith lambda = 1.0; 594a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5955e76c431SBarry Smith lambdaprev = lambda; 5965e76c431SBarry Smith gnormprev = *gnorm; 597329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 598ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 599ddd12767SBarry Smith else lambda = lambdatemp; 600275d25c3SBarry Smith 601e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 60243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 603ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 60443e71028SBarry Smith *flag = PETSC_FALSE; 60543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 60643e71028SBarry Smith goto theend1; 60743e71028SBarry Smith } 608f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6094936397dSBarry Smith if (snes->domainerror) { 6104936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6114936397dSBarry Smith PetscFunctionReturn(0); 61219717074SBarry Smith } 613cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 614f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 61594298653SBarry Smith if (neP->monitor) { 61694298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 61794298653SBarry Smith } 618e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 61994298653SBarry Smith if (neP->monitor) { 62094298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 62194298653SBarry Smith } 622ad922ac9SBarry Smith goto theend1; 6235e76c431SBarry Smith } 6245e76c431SBarry Smith 6255e76c431SBarry Smith /* Fit points with cubic */ 6265e76c431SBarry Smith count = 1; 627bb9195b6SBarry Smith while (PETSC_TRUE) { 628dc357ad2SBarry Smith if (lambda <= minlambda) { 62994298653SBarry Smith if (neP->monitor) { 63094298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 63194298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr); 63294298653SBarry Smith } 63343e71028SBarry Smith *flag = PETSC_FALSE; 63443e71028SBarry Smith break; 6355e76c431SBarry Smith } 6365d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6375d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 638ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6392b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6405e76c431SBarry Smith d = b*b - 3*a*initslope; 6415e76c431SBarry Smith if (d < 0.0) d = 0.0; 6425e76c431SBarry Smith if (a == 0.0) { 6435e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6445e76c431SBarry Smith } else { 6455e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6465e76c431SBarry Smith } 6475e76c431SBarry Smith lambdaprev = lambda; 6485e76c431SBarry Smith gnormprev = *gnorm; 649329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 650329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6515e76c431SBarry Smith else lambda = lambdatemp; 652e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 65343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 654ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 655ae15b995SBarry 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); 656ed98166eSMatthew Knepley *flag = PETSC_FALSE; 65743e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 658ed98166eSMatthew Knepley break; 659ed98166eSMatthew Knepley } 660f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6614936397dSBarry Smith if (snes->domainerror) { 6624936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6634936397dSBarry Smith PetscFunctionReturn(0); 66419717074SBarry Smith } 665cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 666f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 667e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 66894298653SBarry Smith if (neP->monitor) { 66994298653SBarry Smith ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 67094298653SBarry Smith } 67143e71028SBarry Smith break; 6722b022350SLois Curfman McInnes } else { 67394298653SBarry Smith if (neP->monitor) { 67494298653SBarry Smith ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 67594298653SBarry Smith } 6765e76c431SBarry Smith } 6775e76c431SBarry Smith count++; 6785e76c431SBarry Smith } 679ad922ac9SBarry Smith theend1: 6808f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6813c632250SBarry Smith if (neP->postcheckstep && *flag) { 6823c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6833c632250SBarry Smith if (changed_y) { 68479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6853c632250SBarry Smith } 6863c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 687f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6884936397dSBarry Smith if (snes->domainerror) { 6894936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6904936397dSBarry Smith PetscFunctionReturn(0); 69119717074SBarry Smith } 6928f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 693f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 694a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6958f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 696a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6978f99978dSLois Curfman McInnes } 6988f99978dSLois Curfman McInnes } 699d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7003a40ed3dSBarry Smith PetscFunctionReturn(0); 7015e76c431SBarry Smith } 70204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7034a2ae208SSatish Balay #undef __FUNCT__ 70417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 7054b828684SBarry Smith /*@C 70617bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 7075e76c431SBarry Smith 708c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 709c7afd0dbSLois Curfman McInnes 7105e42470aSBarry Smith Input Parameters: 711c7afd0dbSLois Curfman McInnes + snes - the SNES context 712af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 7135e42470aSBarry Smith . x - current iterate 7145e42470aSBarry Smith . f - residual evaluated at x 7153c632250SBarry Smith . y - search direction 7165e42470aSBarry Smith . w - work vector 717dc357ad2SBarry Smith . fnorm - 2-norm of f 718dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 7195e42470aSBarry Smith 720c4a48953SLois Curfman McInnes Output Parameters: 7217f3332b4SBarry Smith + g - residual evaluated at new iterate w 722e106eecfSBarry Smith . w - new iterate (x + lambda*y) 7235e42470aSBarry Smith . gnorm - 2-norm of g 7245e42470aSBarry Smith . ynorm - 2-norm of search length 725a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 726fee21e36SBarry Smith 727ce9499c7SBarry Smith Options Database Keys: 728ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 729ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 730e106eecfSBarry Smith . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 731e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 732e106eecfSBarry Smith 7335e42470aSBarry Smith Notes: 73417bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7355e42470aSBarry Smith 73636851e7fSLois Curfman McInnes Level: advanced 73736851e7fSLois Curfman McInnes 738f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 739f59ffdedSLois Curfman McInnes 74017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7415e42470aSBarry Smith @*/ 742dc357ad2SBarry 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) 7435e76c431SBarry Smith { 744406556e6SLois Curfman McInnes /* 745406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 746406556e6SLois Curfman McInnes minimization problem: 747406556e6SLois Curfman McInnes min z(x): R^n -> R, 748406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 749406556e6SLois Curfman McInnes */ 750e106eecfSBarry Smith PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 751aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 752f5b6597dSBarry Smith PetscScalar cinitslope; 7536b5873e3SBarry Smith #endif 754dfbe8321SBarry Smith PetscErrorCode ierr; 755a7cc72afSBarry Smith PetscInt count; 7564b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 7573c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7585e76c431SBarry Smith 7593a40ed3dSBarry Smith PetscFunctionBegin; 760d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 761a7cc72afSBarry Smith *flag = PETSC_TRUE; 7625e76c431SBarry Smith 7633505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 76463b9fa5eSBarry Smith if (*ynorm == 0.0) { 76594298653SBarry Smith if (neP->monitor) { 76694298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 76794298653SBarry Smith } 768b37302c6SLois Curfman McInnes *gnorm = fnorm; 769e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 770b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 771a7cc72afSBarry Smith *flag = PETSC_FALSE; 772ad922ac9SBarry Smith goto theend2; 77394a9d846SBarry Smith } 774e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 775e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 776e106eecfSBarry Smith *ynorm = neP->maxstep; 7775e76c431SBarry Smith } 778ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 779e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 780a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 781aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 782a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 783329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7845e42470aSBarry Smith #else 785a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7865e42470aSBarry Smith #endif 7875e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7885e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 789f8833479SBarry Smith ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 7905e42470aSBarry Smith 791e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 79243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 793ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 79443e71028SBarry Smith *flag = PETSC_FALSE; 79543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 79643e71028SBarry Smith goto theend2; 79743e71028SBarry Smith } 7984936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 7994936397dSBarry Smith if (snes->domainerror) { 8004936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8014936397dSBarry Smith PetscFunctionReturn(0); 80219717074SBarry Smith } 803cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 804f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 805e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 80694298653SBarry Smith if (neP->monitor) { 80794298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 80894298653SBarry Smith } 809ad922ac9SBarry Smith goto theend2; 8105e42470aSBarry Smith } 8115e42470aSBarry Smith 8125e42470aSBarry Smith /* Fit points with quadratic */ 813313b5538SBarry Smith lambda = 1.0; 8145e42470aSBarry Smith count = 1; 8155ca10a37SBarry Smith while (PETSC_TRUE) { 8165e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 81794298653SBarry Smith if (neP->monitor) { 81894298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 81994298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 82094298653SBarry Smith } 821e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 82243e71028SBarry Smith *flag = PETSC_FALSE; 82343e71028SBarry Smith break; 8245e42470aSBarry Smith } 825a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 826329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 827329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 828329f5518SBarry Smith else lambda = lambdatemp; 829275d25c3SBarry Smith 830e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 83143e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 832ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 833ae15b995SBarry 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); 834ed98166eSMatthew Knepley *flag = PETSC_FALSE; 83543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 836ed98166eSMatthew Knepley break; 837ed98166eSMatthew Knepley } 838f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8394936397dSBarry Smith if (snes->domainerror) { 8404936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8414936397dSBarry Smith PetscFunctionReturn(0); 84219717074SBarry Smith } 843cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 844f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 845e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 84694298653SBarry Smith if (neP->monitor) { 84794298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 84894298653SBarry Smith } 84906259719SBarry Smith break; 8505e42470aSBarry Smith } 8515e42470aSBarry Smith count++; 8525e42470aSBarry Smith } 853ad922ac9SBarry Smith theend2: 8548f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8553c632250SBarry Smith if (neP->postcheckstep) { 8563c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8573c632250SBarry Smith if (changed_y) { 85879f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8593c632250SBarry Smith } 8603c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8613c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 8624936397dSBarry Smith if (snes->domainerror) { 8634936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8644936397dSBarry Smith PetscFunctionReturn(0); 86519717074SBarry Smith } 8668f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 867a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8688f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 869a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 870f66fdb6dSSatish Balay if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8718f99978dSLois Curfman McInnes } 8728f99978dSLois Curfman McInnes } 873d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8743a40ed3dSBarry Smith PetscFunctionReturn(0); 8755e76c431SBarry Smith } 8762343ba6eSBarry Smith 87704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8784a2ae208SSatish Balay #undef __FUNCT__ 87917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 880c9e6a524SLois Curfman McInnes /*@C 88117bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8824b27c08aSLois Curfman McInnes by the method SNESLS. 8835e76c431SBarry Smith 8845e76c431SBarry Smith Input Parameters: 885c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 886af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 887c7afd0dbSLois Curfman McInnes - func - pointer to int function 8885e76c431SBarry Smith 889fee21e36SBarry Smith Collective on SNES 890fee21e36SBarry Smith 891c4a48953SLois Curfman McInnes Available Routines: 89217bae607SBarry Smith + SNESLineSearchCubic() - default line search 89317bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 89417bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 89517bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8965e76c431SBarry Smith 897c4a48953SLois Curfman McInnes Options Database Keys: 8984b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8994b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 900e106eecfSBarry Smith . -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 901e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 902c4a48953SLois Curfman McInnes 9035e76c431SBarry Smith Calling sequence of func: 904bd208895SLois Curfman McInnes .vb 905dc357ad2SBarry 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) 906bd208895SLois Curfman McInnes .ve 9075e76c431SBarry Smith 9085e76c431SBarry Smith Input parameters for func: 909c7afd0dbSLois Curfman McInnes + snes - nonlinear context 910af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 9115e76c431SBarry Smith . x - current iterate 9125e76c431SBarry Smith . f - residual evaluated at x 9133c632250SBarry Smith . y - search direction 914c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9155e76c431SBarry Smith 9165e76c431SBarry Smith Output parameters for func: 917c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9183c632250SBarry Smith . w - new iterate 9195e76c431SBarry Smith . gnorm - 2-norm of g 9205e76c431SBarry Smith . ynorm - 2-norm of search length 921a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 922f59ffdedSLois Curfman McInnes 92336851e7fSLois Curfman McInnes Level: advanced 92436851e7fSLois Curfman McInnes 925f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 926f59ffdedSLois Curfman McInnes 92717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 92817bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9295e76c431SBarry Smith @*/ 930dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 9315e76c431SBarry Smith { 932dc357ad2SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 93382bf6240SBarry Smith 9343a40ed3dSBarry Smith PetscFunctionBegin; 93517bae607SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 93682bf6240SBarry Smith if (f) { 937af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 93882bf6240SBarry Smith } 9393a40ed3dSBarry Smith PetscFunctionReturn(0); 9405e76c431SBarry Smith } 9418e019c35SBarry Smith 942dc357ad2SBarry 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*/ 94304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 944fb2e594dSBarry Smith EXTERN_C_BEGIN 9454a2ae208SSatish Balay #undef __FUNCT__ 94617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 94717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 94882bf6240SBarry Smith { 94982bf6240SBarry Smith PetscFunctionBegin; 9504b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9514b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 95282bf6240SBarry Smith PetscFunctionReturn(0); 95382bf6240SBarry Smith } 954fb2e594dSBarry Smith EXTERN_C_END 95504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9564a2ae208SSatish Balay #undef __FUNCT__ 9573c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 958c8dd0c0dSLois Curfman McInnes /*@C 9593c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9604b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 961c8dd0c0dSLois Curfman McInnes 962c8dd0c0dSLois Curfman McInnes Input Parameters: 963c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9643c632250SBarry Smith . func - pointer to function 965c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 966c8dd0c0dSLois Curfman McInnes 967c8dd0c0dSLois Curfman McInnes Collective on SNES 968c8dd0c0dSLois Curfman McInnes 969c8dd0c0dSLois Curfman McInnes Calling sequence of func: 970c8dd0c0dSLois Curfman McInnes .vb 9713c632250SBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 972c8dd0c0dSLois Curfman McInnes .ve 973b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 974b1ae27eaSLois Curfman McInnes on failure. 975c8dd0c0dSLois Curfman McInnes 976c8dd0c0dSLois Curfman McInnes Input parameters for func: 977c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 978c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9793c632250SBarry Smith . x - previous iterate 9803c632250SBarry Smith . y - new search direction and length 9813c632250SBarry Smith - w - current candidate iterate 982c8dd0c0dSLois Curfman McInnes 983c8dd0c0dSLois Curfman McInnes Output parameters for func: 9843c632250SBarry Smith + y - search direction (possibly changed) 9853c632250SBarry Smith . w - current iterate (possibly modified) 9863c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 9873c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 988c8dd0c0dSLois Curfman McInnes 989c8dd0c0dSLois Curfman McInnes Level: advanced 990c8dd0c0dSLois Curfman McInnes 9919e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 9929e247f21SBarry Smith 9933c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 9943c632250SBarry Smith 9953c632250SBarry Smith On input w = x + y 9963c632250SBarry Smith 99717bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 998b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 999ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 10008f99978dSLois Curfman McInnes 100117bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 1002ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 1003ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 1004ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 10059e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 1006b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 10078f99978dSLois Curfman McInnes 1008c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 1009c8dd0c0dSLois Curfman McInnes 101017bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1011c8dd0c0dSLois Curfman McInnes @*/ 10123c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 1013c8dd0c0dSLois Curfman McInnes { 10143c632250SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 1015c8dd0c0dSLois Curfman McInnes 1016c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10173c632250SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1018c8dd0c0dSLois Curfman McInnes if (f) { 1019c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1020c8dd0c0dSLois Curfman McInnes } 1021c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1022c8dd0c0dSLois Curfman McInnes } 102394298653SBarry Smith 10249c3ca13aSBarry Smith #undef __FUNCT__ 10259c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10269c3ca13aSBarry Smith /*@C 10279c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10287e4bb74cSBarry Smith before the line search is called. 10299c3ca13aSBarry Smith 10309c3ca13aSBarry Smith Input Parameters: 10319c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10329c3ca13aSBarry Smith . func - pointer to function 10339c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10349c3ca13aSBarry Smith 10359c3ca13aSBarry Smith Collective on SNES 10369c3ca13aSBarry Smith 10379c3ca13aSBarry Smith Calling sequence of func: 10389c3ca13aSBarry Smith .vb 10391e3347e8SBarry Smith int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 10409c3ca13aSBarry Smith .ve 10419c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10429c3ca13aSBarry Smith on failure. 10439c3ca13aSBarry Smith 10449c3ca13aSBarry Smith Input parameters for func: 10459c3ca13aSBarry Smith + snes - nonlinear context 10469c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10479c3ca13aSBarry Smith . x - previous iterate 10489c3ca13aSBarry Smith - y - new search direction and length 10499c3ca13aSBarry Smith 10509c3ca13aSBarry Smith Output parameters for func: 10519c3ca13aSBarry Smith + y - search direction (possibly changed) 10529c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10539c3ca13aSBarry Smith 10549c3ca13aSBarry Smith Level: advanced 10559c3ca13aSBarry Smith 10569c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10579c3ca13aSBarry Smith 10589c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10599c3ca13aSBarry Smith 10607e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 10619c3ca13aSBarry Smith @*/ 10629c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 10639c3ca13aSBarry Smith { 10649c3ca13aSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 10659c3ca13aSBarry Smith 10669c3ca13aSBarry Smith PetscFunctionBegin; 10679c3ca13aSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 10689c3ca13aSBarry Smith if (f) { 10699c3ca13aSBarry Smith ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 10709c3ca13aSBarry Smith } 10719c3ca13aSBarry Smith PetscFunctionReturn(0); 10729c3ca13aSBarry Smith } 10739c3ca13aSBarry Smith 107494298653SBarry Smith #undef __FUNCT__ 107594298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor" 107694298653SBarry Smith /*@C 107794298653SBarry Smith SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search 107894298653SBarry Smith 107994298653SBarry Smith Input Parameters: 108094298653SBarry Smith + snes - nonlinear context obtained from SNESCreate() 108194298653SBarry Smith - flg - PETSC_TRUE to monitor the line search 108294298653SBarry Smith 108394298653SBarry Smith Collective on SNES 108494298653SBarry Smith 108594298653SBarry Smith Options Database: 108694298653SBarry Smith . -snes_ls_monitor 108794298653SBarry Smith 108894298653SBarry Smith Level: intermediate 108994298653SBarry Smith 109094298653SBarry Smith 109194298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 109294298653SBarry Smith @*/ 109394298653SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor(SNES snes,PetscTruth flg) 109494298653SBarry Smith { 109594298653SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscTruth); 109694298653SBarry Smith 109794298653SBarry Smith PetscFunctionBegin; 109894298653SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetMonitor_C",(void (**)(void))&f);CHKERRQ(ierr); 109994298653SBarry Smith if (f) { 110094298653SBarry Smith ierr = (*f)(snes,flg);CHKERRQ(ierr); 110194298653SBarry Smith } 110294298653SBarry Smith PetscFunctionReturn(0); 110394298653SBarry Smith } 110494298653SBarry Smith 1105c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 11069c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 11079c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1108c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 11094a2ae208SSatish Balay #undef __FUNCT__ 11103c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 11119c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1112c8dd0c0dSLois Curfman McInnes { 1113c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 11143c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 11153c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1116c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1117c8dd0c0dSLois Curfman McInnes } 1118c8dd0c0dSLois Curfman McInnes EXTERN_C_END 11199c3ca13aSBarry Smith 11209c3ca13aSBarry Smith EXTERN_C_BEGIN 11219c3ca13aSBarry Smith #undef __FUNCT__ 11229c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 11239c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 11249c3ca13aSBarry Smith { 11259c3ca13aSBarry Smith PetscFunctionBegin; 11269c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 11279c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 11289c3ca13aSBarry Smith PetscFunctionReturn(0); 11299c3ca13aSBarry Smith } 11309c3ca13aSBarry Smith EXTERN_C_END 1131329e7e40SMatthew Knepley 113294298653SBarry Smith EXTERN_C_BEGIN 113394298653SBarry Smith #undef __FUNCT__ 113494298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS" 113594298653SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_LS(SNES snes,PetscTruth flg) 113694298653SBarry Smith { 113794298653SBarry Smith SNES_LS *ls = (SNES_LS*)snes->data; 113894298653SBarry Smith PetscErrorCode ierr; 113994298653SBarry Smith 114094298653SBarry Smith PetscFunctionBegin; 114194298653SBarry Smith if (flg && !ls->monitor) { 114294298653SBarry Smith ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&ls->monitor);CHKERRQ(ierr); 114394298653SBarry Smith } else if (!flg && ls->monitor) { 114494298653SBarry Smith ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 114594298653SBarry Smith } 114694298653SBarry Smith PetscFunctionReturn(0); 114794298653SBarry Smith } 114894298653SBarry Smith EXTERN_C_END 114994298653SBarry Smith 1150329e7e40SMatthew Knepley /* 11514b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 115204d965bbSLois Curfman McInnes 115304d965bbSLois Curfman McInnes Input Parameters: 115404d965bbSLois Curfman McInnes . SNES - the SNES context 115504d965bbSLois Curfman McInnes . viewer - visualization context 115604d965bbSLois Curfman McInnes 115704d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 115804d965bbSLois Curfman McInnes */ 11594a2ae208SSatish Balay #undef __FUNCT__ 11604b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11616849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1162a935fc98SLois Curfman McInnes { 11634b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11642fc52814SBarry Smith const char *cstr; 1165dfbe8321SBarry Smith PetscErrorCode ierr; 116632077d6dSBarry Smith PetscTruth iascii; 1167a935fc98SLois Curfman McInnes 11683a40ed3dSBarry Smith PetscFunctionBegin; 116932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 117032077d6dSBarry Smith if (iascii) { 117117bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 117217bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 117317bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 117419bcc07fSBarry Smith else cstr = "unknown"; 1175b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1176e106eecfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 11775cd90555SBarry Smith } else { 117879a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 117919bcc07fSBarry Smith } 11803a40ed3dSBarry Smith PetscFunctionReturn(0); 1181a935fc98SLois Curfman McInnes } 118204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 118304d965bbSLois Curfman McInnes /* 11844b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 118504d965bbSLois Curfman McInnes 118604d965bbSLois Curfman McInnes Input Parameter: 118704d965bbSLois Curfman McInnes . snes - the SNES context 118804d965bbSLois Curfman McInnes 118904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 119004d965bbSLois Curfman McInnes */ 11914a2ae208SSatish Balay #undef __FUNCT__ 11924b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11936849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11945e42470aSBarry Smith { 11954b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1196e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1197dfbe8321SBarry Smith PetscErrorCode ierr; 1198a7cc72afSBarry Smith PetscInt indx; 119994298653SBarry Smith PetscTruth flg,set; 12005e42470aSBarry Smith 12013a40ed3dSBarry Smith PetscFunctionBegin; 1202b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 12034b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 12044b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1205e106eecfSBarry Smith ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 120694298653SBarry Smith ierr = PetscOptionsTruth("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 120794298653SBarry Smith if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1208186905e3SBarry Smith 120917bae607SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 121025cdf11fSBarry Smith if (flg) { 121122e36657SBarry Smith switch (indx) { 1212b49fd9e1SBarry Smith case 0: 121317bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1214b49fd9e1SBarry Smith break; 1215b49fd9e1SBarry Smith case 1: 121617bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1217b49fd9e1SBarry Smith break; 1218b49fd9e1SBarry Smith case 2: 121917bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1220b49fd9e1SBarry Smith break; 1221b49fd9e1SBarry Smith case 3: 122217bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1223b49fd9e1SBarry Smith break; 12245e42470aSBarry Smith } 12255e42470aSBarry Smith } 1226b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 12273a40ed3dSBarry Smith PetscFunctionReturn(0); 12285e42470aSBarry Smith } 122904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1230ebe3b25bSBarry Smith /*MC 1231ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 123204d965bbSLois Curfman McInnes 1233ebe3b25bSBarry Smith Options Database: 1234ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1235ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1236e106eecfSBarry Smith . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 1237*acbee50cSBarry Smith . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1238*acbee50cSBarry Smith - -snes_ls_monitor - print information about progress of line searches 1239*acbee50cSBarry Smith 124004d965bbSLois Curfman McInnes 1241ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 124204d965bbSLois Curfman McInnes 1243ee3001cbSBarry Smith Level: beginner 1244ee3001cbSBarry Smith 124517bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 124617bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 124717bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1248ebe3b25bSBarry Smith 1249ebe3b25bSBarry Smith M*/ 1250fb2e594dSBarry Smith EXTERN_C_BEGIN 12514a2ae208SSatish Balay #undef __FUNCT__ 12524b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 125363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 12545e42470aSBarry Smith { 1255dfbe8321SBarry Smith PetscErrorCode ierr; 12564b27c08aSLois Curfman McInnes SNES_LS *neP; 12575e42470aSBarry Smith 12583a40ed3dSBarry Smith PetscFunctionBegin; 1259e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 1260e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 1261e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 1262e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 1263e7788613SBarry Smith snes->ops->view = SNESView_LS; 12645e42470aSBarry Smith 126538f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 12665e42470aSBarry Smith snes->data = (void*)neP; 12675e42470aSBarry Smith neP->alpha = 1.e-4; 12685e42470aSBarry Smith neP->maxstep = 1.e8; 1269e106eecfSBarry Smith neP->minlambda = 1.e-12; 127017bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1271c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12723c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12733c632250SBarry Smith neP->postcheck = PETSC_NULL; 12743c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12753c632250SBarry Smith neP->precheck = PETSC_NULL; 127682bf6240SBarry Smith 127794298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr); 127894298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 127994298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 128094298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 128182bf6240SBarry Smith 12823a40ed3dSBarry Smith PetscFunctionReturn(0); 12835e42470aSBarry Smith } 1284fb2e594dSBarry Smith EXTERN_C_END 128582bf6240SBarry Smith 128682bf6240SBarry Smith 128782bf6240SBarry Smith 1288