15e76c431SBarry Smith 2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h> 35e42470aSBarry Smith 48a5d9ee4SBarry Smith /* 520f52e01SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 620f52e01SBarry 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 736669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 820f52e01SBarry Smith for this trick. One assumes that the probability that W is in the null space of J is very, very small. 98a5d9ee4SBarry Smith */ 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 12ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 138a5d9ee4SBarry Smith { 148a5d9ee4SBarry Smith PetscReal a1; 15dfbe8321SBarry Smith PetscErrorCode ierr; 16ace3abfcSBarry Smith PetscBool hastranspose; 178a5d9ee4SBarry Smith 188a5d9ee4SBarry Smith PetscFunctionBegin; 198a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2036669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2136669109SBarry Smith if (hastranspose) { 228a5d9ee4SBarry Smith /* Compute || J^T F|| */ 238a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 248a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 251e2582c4SBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 268a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2736669109SBarry Smith } else { 2836669109SBarry Smith Vec work; 29ea709b57SSatish Balay PetscScalar result; 3036669109SBarry Smith PetscReal wnorm; 3136669109SBarry Smith 32abb0e124SMatthew Knepley ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 3336669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3536669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 376bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 3836669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 391e2582c4SBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 4036669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4136669109SBarry Smith } 428a5d9ee4SBarry Smith PetscFunctionReturn(0); 438a5d9ee4SBarry Smith } 448a5d9ee4SBarry Smith 4574637425SBarry Smith /* 465ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4774637425SBarry Smith */ 484a2ae208SSatish Balay #undef __FUNCT__ 494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 501e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 5174637425SBarry Smith { 5274637425SBarry Smith PetscReal a1,a2; 53dfbe8321SBarry Smith PetscErrorCode ierr; 54ace3abfcSBarry Smith PetscBool hastranspose; 5574637425SBarry Smith 5674637425SBarry Smith PetscFunctionBegin; 5774637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 5874637425SBarry Smith if (hastranspose) { 5974637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6079f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6174637425SBarry Smith 6274637425SBarry Smith /* Compute || J^T W|| */ 6374637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6474637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6574637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 6675567043SBarry Smith if (a1 != 0.0) { 671e2582c4SBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 6885471664SBarry Smith } 6974637425SBarry Smith } 7074637425SBarry Smith PetscFunctionReturn(0); 7174637425SBarry Smith } 7274637425SBarry Smith 7304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7404d965bbSLois Curfman McInnes 7504d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7694b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7704d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7804d965bbSLois Curfman McInnes respectively. 7904d965bbSLois Curfman McInnes 80fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8104d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8204d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 83fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8404d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8504d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 864b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 874b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 8804d965bbSLois Curfman McInnes These routines are actually called via the common user interface 8904d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9004d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9104d965bbSLois Curfman McInnes for all nonlinear solvers. 9204d965bbSLois Curfman McInnes 9304d965bbSLois Curfman McInnes Another key routine is: 9404d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9504d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9604d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9704d965bbSLois Curfman McInnes SNESSolve() if necessary. 9804d965bbSLois Curfman McInnes 9904d965bbSLois Curfman McInnes Additional basic routines are: 10004d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10104d965bbSLois Curfman McInnes have actually been used. 10204d965bbSLois Curfman McInnes These are called by application codes via the interface routines 103186905e3SBarry Smith SNESView(). 10404d965bbSLois Curfman McInnes 10504d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10604d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10704d965bbSLois Curfman McInnes above description applies to these categories also. 10804d965bbSLois Curfman McInnes 10904d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1105e42470aSBarry Smith /* 1114b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11204d965bbSLois Curfman McInnes method with a line search. 1135e76c431SBarry Smith 11404d965bbSLois Curfman McInnes Input Parameters: 11504d965bbSLois Curfman McInnes . snes - the SNES context 1165e76c431SBarry Smith 11704d965bbSLois Curfman McInnes Output Parameter: 118c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 11904d965bbSLois Curfman McInnes 12004d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1215e76c431SBarry Smith 1225e76c431SBarry Smith Notes: 1235e76c431SBarry Smith This implements essentially a truncated Newton method with a 1245e76c431SBarry Smith line search. By default a cubic backtracking line search 1255e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1265e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 127393d2d9aSLois Curfman McInnes and Schnabel. 1285e76c431SBarry Smith */ 1294a2ae208SSatish Balay #undef __FUNCT__ 1304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 131dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1325e76c431SBarry Smith { 1334b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1346849ba73SBarry Smith PetscErrorCode ierr; 135a7cc72afSBarry Smith PetscInt maxits,i,lits; 136ace3abfcSBarry Smith PetscBool lssucceed; 137112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 13885385478SLisandro Dalcin PetscReal fnorm,gnorm,xnorm=0,ynorm; 13985385478SLisandro Dalcin Vec Y,X,F,G,W; 1403d4c4710SBarry Smith KSPConvergedReason kspreason; 1415e76c431SBarry Smith 1423a40ed3dSBarry Smith PetscFunctionBegin; 1433d4c4710SBarry Smith snes->numFailures = 0; 1443d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 145184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 146184914b5SBarry Smith 1475e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1485e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1505e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1515e42470aSBarry Smith G = snes->work[1]; 1525e42470aSBarry Smith W = snes->work[2]; 1535e76c431SBarry Smith 1544c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1554c49b128SBarry Smith snes->iter = 0; 15675567043SBarry Smith snes->norm = 0.0; 1574c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15819717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1594936397dSBarry Smith if (snes->domainerror) { 1604936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1614936397dSBarry Smith PetscFunctionReturn(0); 1624936397dSBarry Smith } 1632613ca53SBarry Smith ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1642613ca53SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1652613ca53SBarry Smith ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 1662613ca53SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 16762cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1680f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1695e42470aSBarry Smith snes->norm = fnorm; 1700f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 17184c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 172*7a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1733f149594SLisandro Dalcin 1743f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1753f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 17685385478SLisandro Dalcin /* test convergence */ 17785385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 17806ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 179d034289bSBarry Smith 1805e76c431SBarry Smith for (i=0; i<maxits; i++) { 1815e76c431SBarry Smith 182329e7e40SMatthew Knepley /* Call general purpose update function */ 183e7788613SBarry Smith if (snes->ops->update) { 184e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 185de8cb200SMatthew Knepley } 186329e7e40SMatthew Knepley 187ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1885334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 18994b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 19071f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 1913d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1923d4c4710SBarry Smith if (kspreason < 0) { 1933d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 194ef998cc9SBarry 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); 1953d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 196ab81109fSSatish Balay break; 1973d4c4710SBarry Smith } 1983d4c4710SBarry Smith } 1993d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2003f149594SLisandro Dalcin snes->linear_its += lits; 2013f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 20274637425SBarry Smith 2039c3ca13aSBarry Smith if (neP->precheckstep) { 204ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 2059c3ca13aSBarry Smith ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 2069c3ca13aSBarry Smith } 2079c3ca13aSBarry Smith 208b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 2091e2582c4SBarry Smith ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 21085471664SBarry Smith } 21174637425SBarry Smith 212ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 213ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 214e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 215ea4d3ed3SLois Curfman McInnes */ 21685385478SLisandro Dalcin ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 2173f149594SLisandro Dalcin ynorm = 1; gnorm = fnorm; 218dc357ad2SBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 219ae15b995SBarry Smith ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 2204a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2214936397dSBarry Smith if (snes->domainerror) { 2224936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2234936397dSBarry Smith PetscFunctionReturn(0); 2244936397dSBarry Smith } 225a7cc72afSBarry Smith if (!lssucceed) { 22650ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 227ace3abfcSBarry Smith PetscBool ismin; 228647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2291e2582c4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 2308a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2313505fcc1SBarry Smith break; 2323505fcc1SBarry Smith } 23350ffb88aSMatthew Knepley } 23485385478SLisandro Dalcin /* Update function and solution vectors */ 23585385478SLisandro Dalcin fnorm = gnorm; 23685385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 23785385478SLisandro Dalcin ierr = VecCopy(W,X);CHKERRQ(ierr); 23885385478SLisandro Dalcin /* Monitor convergence */ 23985385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 24085385478SLisandro Dalcin snes->iter = i+1; 24185385478SLisandro Dalcin snes->norm = fnorm; 24285385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24385385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 244*7a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 24585385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 24685385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 247e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2483f149594SLisandro Dalcin if (snes->reason) break; 24916c95cacSBarry Smith } 25052392280SLois Curfman McInnes if (i == maxits) { 251ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25285385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25352392280SLois Curfman McInnes } 2543a40ed3dSBarry Smith PetscFunctionReturn(0); 2555e76c431SBarry Smith } 25604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25704d965bbSLois Curfman McInnes /* 2584b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2594b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26004d965bbSLois Curfman McInnes 26104d965bbSLois Curfman McInnes Input Parameter: 26204d965bbSLois Curfman McInnes . snes - the SNES context 26304d965bbSLois Curfman McInnes . x - the solution vector 26404d965bbSLois Curfman McInnes 26504d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26604d965bbSLois Curfman McInnes 26704d965bbSLois Curfman McInnes Notes: 2684b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 26904d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27004d965bbSLois Curfman McInnes the call to SNESSolve(). 27104d965bbSLois Curfman McInnes */ 2724a2ae208SSatish Balay #undef __FUNCT__ 2734b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 274dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2755e76c431SBarry Smith { 276dfbe8321SBarry Smith PetscErrorCode ierr; 2773a40ed3dSBarry Smith 2783a40ed3dSBarry Smith PetscFunctionBegin; 27958c9b817SLisandro Dalcin ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 2803a40ed3dSBarry Smith PetscFunctionReturn(0); 2815e76c431SBarry Smith } 28204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2836b8b9a38SLisandro Dalcin 2846b8b9a38SLisandro Dalcin #undef __FUNCT__ 2856b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS" 2866b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes) 2876b8b9a38SLisandro Dalcin { 2886b8b9a38SLisandro Dalcin PetscErrorCode ierr; 2896b8b9a38SLisandro Dalcin 2906b8b9a38SLisandro Dalcin PetscFunctionBegin; 2916b8b9a38SLisandro Dalcin if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2926b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2936b8b9a38SLisandro Dalcin } 2946b8b9a38SLisandro Dalcin 29504d965bbSLois Curfman McInnes /* 2964b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2974b27c08aSLois Curfman McInnes with SNESCreate_LS(). 29804d965bbSLois Curfman McInnes 29904d965bbSLois Curfman McInnes Input Parameter: 30004d965bbSLois Curfman McInnes . snes - the SNES context 30104d965bbSLois Curfman McInnes 302de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30304d965bbSLois Curfman McInnes */ 3044a2ae208SSatish Balay #undef __FUNCT__ 3054b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 306dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3075e76c431SBarry Smith { 30894298653SBarry Smith SNES_LS *ls = (SNES_LS*) snes->data; 309dfbe8321SBarry Smith PetscErrorCode ierr; 3103a40ed3dSBarry Smith 3113a40ed3dSBarry Smith PetscFunctionBegin; 3126b8b9a38SLisandro Dalcin ierr = SNESReset_LS(snes);CHKERRQ(ierr); 3136bf464f9SBarry Smith ierr = PetscViewerASCIIMonitorDestroy(&ls->monitor);CHKERRQ(ierr); 3145c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 315557d3f75SLisandro Dalcin 316557d3f75SLisandro Dalcin /* clear composed functions */ 31758c9b817SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 31858c9b817SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);CHKERRQ(ierr); 319557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 320557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 321557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 322557d3f75SLisandro Dalcin 3233a40ed3dSBarry Smith PetscFunctionReturn(0); 3245e76c431SBarry Smith } 32504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3264a2ae208SSatish Balay #undef __FUNCT__ 32717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 32882bf6240SBarry Smith 3294b828684SBarry Smith /*@C 33017bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3315e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3325e76c431SBarry Smith to serve as a template and is not recommended for general use. 3335e76c431SBarry Smith 3343f9fe445SBarry Smith Logically Collective on SNES and Vec 335c7afd0dbSLois Curfman McInnes 3365e76c431SBarry Smith Input Parameters: 337c7afd0dbSLois Curfman McInnes + snes - nonlinear context 338af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3395e76c431SBarry Smith . x - current iterate 3405e76c431SBarry Smith . f - residual evaluated at x 3413c632250SBarry Smith . y - search direction 342dc357ad2SBarry Smith . fnorm - 2-norm of f 343dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 3445e76c431SBarry Smith 345c4a48953SLois Curfman McInnes Output Parameters: 346c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3473c632250SBarry Smith . w - new iterate 3485e76c431SBarry Smith . gnorm - 2-norm of g 3495e76c431SBarry Smith . ynorm - 2-norm of search length 350a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 351fee21e36SBarry Smith 352c4a48953SLois Curfman McInnes Options Database Key: 35317bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 354c4a48953SLois Curfman McInnes 35536851e7fSLois Curfman McInnes Level: advanced 35636851e7fSLois Curfman McInnes 35728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 35828ae5a21SLois Curfman McInnes 35917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 36017bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3615e76c431SBarry Smith @*/ 3627087cfbeSBarry Smith PetscErrorCode SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 3635e76c431SBarry Smith { 364dfbe8321SBarry Smith PetscErrorCode ierr; 3654b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 366ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3675334005bSBarry Smith 3683a40ed3dSBarry Smith PetscFunctionBegin; 369a7cc72afSBarry Smith *flag = PETSC_TRUE; 370d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 371a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 37279f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3733c632250SBarry Smith if (neP->postcheckstep) { 3743c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3758f99978dSLois Curfman McInnes } 3763c632250SBarry Smith if (changed_y) { 37779f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3783c632250SBarry Smith } 3794936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 3804936397dSBarry Smith if (!snes->domainerror) { 381a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 38262cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3834936397dSBarry Smith } 384d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3853a40ed3dSBarry Smith PetscFunctionReturn(0); 3865e76c431SBarry Smith } 38704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3884a2ae208SSatish Balay #undef __FUNCT__ 38917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 39082bf6240SBarry Smith 39129e0b56fSBarry Smith /*@C 39217bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 39329e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 39429e0b56fSBarry Smith even compute the norm of the function or search direction; this 39529e0b56fSBarry Smith is intended only when you know the full step is fine and are 39629e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 397c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 398c7afd0dbSLois Curfman McInnes 3993f9fe445SBarry Smith Logically Collective on SNES and Vec 40029e0b56fSBarry Smith 40129e0b56fSBarry Smith Input Parameters: 402c7afd0dbSLois Curfman McInnes + snes - nonlinear context 403af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 40429e0b56fSBarry Smith . x - current iterate 40529e0b56fSBarry Smith . f - residual evaluated at x 4063c632250SBarry Smith . y - search direction 40729e0b56fSBarry Smith . w - work vector 408dc357ad2SBarry Smith . fnorm - 2-norm of f 409dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 41029e0b56fSBarry Smith 41129e0b56fSBarry Smith Output Parameters: 412c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4133c632250SBarry Smith . w - new iterate 41429e0b56fSBarry Smith . gnorm - not changed 41529e0b56fSBarry Smith . ynorm - not changed 416a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 417fee21e36SBarry Smith 41829e0b56fSBarry Smith Options Database Key: 41917bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 42029e0b56fSBarry Smith 4218cbba510SBarry Smith Notes: 42217bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 423ea56f5baSLois Curfman McInnes either the options 424ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 425ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 426329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 427329f5518SBarry Smith otherwise, the SNES solver will generate an error. 428329f5518SBarry Smith 429329f5518SBarry Smith During the final iteration this will not evaluate the function at 430329f5518SBarry Smith the solution point. This is to save a function evaluation while 431329f5518SBarry Smith using pseudo-timestepping. 4328cbba510SBarry Smith 433ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 434a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 435ea56f5baSLois Curfman McInnes correct, since they are not computed. 436ea56f5baSLois Curfman McInnes 437ea56f5baSLois Curfman McInnes Level: advanced 4388cbba510SBarry Smith 43929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 44029e0b56fSBarry Smith 44117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 44217bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 44329e0b56fSBarry Smith @*/ 4447087cfbeSBarry Smith PetscErrorCode SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 44529e0b56fSBarry Smith { 446dfbe8321SBarry Smith PetscErrorCode ierr; 4474b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 448ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 44929e0b56fSBarry Smith 4503a40ed3dSBarry Smith PetscFunctionBegin; 451a7cc72afSBarry Smith *flag = PETSC_TRUE; 452d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 45379f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4543c632250SBarry Smith if (neP->postcheckstep) { 4553c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4563c632250SBarry Smith } 4573c632250SBarry Smith if (changed_y) { 45879f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4598f99978dSLois Curfman McInnes } 460329f5518SBarry Smith 461329f5518SBarry Smith /* don't evaluate function the last time through */ 462329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4634936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 464329f5518SBarry Smith } 465d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4663a40ed3dSBarry Smith PetscFunctionReturn(0); 46729e0b56fSBarry Smith } 46804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4694a2ae208SSatish Balay #undef __FUNCT__ 47017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4714b828684SBarry Smith /*@C 47217bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4735e76c431SBarry Smith 474c7afd0dbSLois Curfman McInnes Collective on SNES 475c7afd0dbSLois Curfman McInnes 4765e76c431SBarry Smith Input Parameters: 477c7afd0dbSLois Curfman McInnes + snes - nonlinear context 478af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4795e76c431SBarry Smith . x - current iterate 4805e76c431SBarry Smith . f - residual evaluated at x 4813c632250SBarry Smith . y - search direction 4825e76c431SBarry Smith . w - work vector 483dc357ad2SBarry Smith . fnorm - 2-norm of f 484dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 4855e76c431SBarry Smith 486393d2d9aSLois Curfman McInnes Output Parameters: 487c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4883c632250SBarry Smith . w - new iterate 4895e76c431SBarry Smith . gnorm - 2-norm of g 4905e76c431SBarry Smith . ynorm - 2-norm of search length 491a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 492fee21e36SBarry Smith 493c4a48953SLois Curfman McInnes Options Database Key: 494e106eecfSBarry Smith + -snes_ls cubic - Activates SNESLineSearchCubic() 495e106eecfSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 496e106eecfSBarry 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) 497e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 498e106eecfSBarry Smith 499c4a48953SLois Curfman McInnes 5005e76c431SBarry Smith Notes: 5015e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 5025e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 5035e76c431SBarry Smith 50436851e7fSLois Curfman McInnes Level: advanced 50536851e7fSLois Curfman McInnes 50628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 50728ae5a21SLois Curfman McInnes 50817bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 5095e76c431SBarry Smith @*/ 5107087cfbeSBarry Smith PetscErrorCode SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 5115e76c431SBarry Smith { 512406556e6SLois Curfman McInnes /* 513406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 514406556e6SLois Curfman McInnes minimization problem: 515406556e6SLois Curfman McInnes min z(x): R^n -> R, 516406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 517406556e6SLois Curfman McInnes */ 518406556e6SLois Curfman McInnes 519e106eecfSBarry Smith PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 520e106eecfSBarry Smith PetscReal minlambda,lambda,lambdatemp; 521aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 522f5b6597dSBarry Smith PetscScalar cinitslope; 5236b5873e3SBarry Smith #endif 524dfbe8321SBarry Smith PetscErrorCode ierr; 525a7cc72afSBarry Smith PetscInt count; 5264b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 527ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 52894298653SBarry Smith MPI_Comm comm; 5295e76c431SBarry Smith 5303a40ed3dSBarry Smith PetscFunctionBegin; 53194298653SBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 532d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 533a7cc72afSBarry Smith *flag = PETSC_TRUE; 5345e76c431SBarry Smith 535cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 53675567043SBarry Smith if (*ynorm == 0.0) { 53794298653SBarry Smith if (neP->monitor) { 53894298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 53994298653SBarry Smith } 540a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5413c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 542a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 543a7cc72afSBarry Smith *flag = PETSC_FALSE; 544ad922ac9SBarry Smith goto theend1; 54594a9d846SBarry Smith } 546e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 54794298653SBarry Smith if (neP->monitor) { 54894298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 54994298653SBarry Smith } 550e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 551e106eecfSBarry Smith *ynorm = neP->maxstep; 5525e76c431SBarry Smith } 553ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 554e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 555a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 556aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 557a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 558329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5595e42470aSBarry Smith #else 560a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5615e42470aSBarry Smith #endif 5625e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5635e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5645e76c431SBarry Smith 565e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 56643e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 567ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 56843e71028SBarry Smith *flag = PETSC_FALSE; 56943e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 57043e71028SBarry Smith goto theend1; 57143e71028SBarry Smith } 5724936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5734936397dSBarry Smith if (snes->domainerror) { 5744936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5754936397dSBarry Smith PetscFunctionReturn(0); 57619717074SBarry Smith } 577313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 57862cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 579f5b6597dSBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 580e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 58194298653SBarry Smith if (neP->monitor) { 58294298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 58394298653SBarry Smith } 584ad922ac9SBarry Smith goto theend1; 5855e76c431SBarry Smith } 5865e76c431SBarry Smith 5875e76c431SBarry Smith /* Fit points with quadratic */ 588313b5538SBarry Smith lambda = 1.0; 589a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5905e76c431SBarry Smith lambdaprev = lambda; 5915e76c431SBarry Smith gnormprev = *gnorm; 592329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 593ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 594ddd12767SBarry Smith else lambda = lambdatemp; 595275d25c3SBarry Smith 596e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 59743e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 598ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 59943e71028SBarry Smith *flag = PETSC_FALSE; 60043e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 60143e71028SBarry Smith goto theend1; 60243e71028SBarry Smith } 603f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6044936397dSBarry Smith if (snes->domainerror) { 6054936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6064936397dSBarry Smith PetscFunctionReturn(0); 60719717074SBarry Smith } 608cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 60962cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 61094298653SBarry Smith if (neP->monitor) { 61194298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 61294298653SBarry Smith } 613e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 61494298653SBarry Smith if (neP->monitor) { 61594298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 61694298653SBarry Smith } 617ad922ac9SBarry Smith goto theend1; 6185e76c431SBarry Smith } 6195e76c431SBarry Smith 6205e76c431SBarry Smith /* Fit points with cubic */ 6215e76c431SBarry Smith count = 1; 622bb9195b6SBarry Smith while (PETSC_TRUE) { 623dc357ad2SBarry Smith if (lambda <= minlambda) { 62494298653SBarry Smith if (neP->monitor) { 62594298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 626d9822059SBarry 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",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr); 62794298653SBarry Smith } 62843e71028SBarry Smith *flag = PETSC_FALSE; 62943e71028SBarry Smith break; 6305e76c431SBarry Smith } 6315d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6325d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 633ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6342b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6355e76c431SBarry Smith d = b*b - 3*a*initslope; 6365e76c431SBarry Smith if (d < 0.0) d = 0.0; 6375e76c431SBarry Smith if (a == 0.0) { 6385e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6395e76c431SBarry Smith } else { 6405e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6415e76c431SBarry Smith } 6425e76c431SBarry Smith lambdaprev = lambda; 6435e76c431SBarry Smith gnormprev = *gnorm; 644329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 645329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6465e76c431SBarry Smith else lambda = lambdatemp; 647e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 64843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 649ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 650ae15b995SBarry 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); 651ed98166eSMatthew Knepley *flag = PETSC_FALSE; 65243e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 653ed98166eSMatthew Knepley break; 654ed98166eSMatthew Knepley } 655f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6564936397dSBarry Smith if (snes->domainerror) { 6574936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6584936397dSBarry Smith PetscFunctionReturn(0); 65919717074SBarry Smith } 660cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 66162cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 662e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 66394298653SBarry Smith if (neP->monitor) { 664d9822059SBarry Smith ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr); 66594298653SBarry Smith } 66643e71028SBarry Smith break; 6672b022350SLois Curfman McInnes } else { 66894298653SBarry Smith if (neP->monitor) { 669d9822059SBarry Smith ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr); 67094298653SBarry Smith } 6715e76c431SBarry Smith } 6725e76c431SBarry Smith count++; 6735e76c431SBarry Smith } 674ad922ac9SBarry Smith theend1: 6758f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6763c632250SBarry Smith if (neP->postcheckstep && *flag) { 6773c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6783c632250SBarry Smith if (changed_y) { 67979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6803c632250SBarry Smith } 6813c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 682f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6834936397dSBarry Smith if (snes->domainerror) { 6844936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6854936397dSBarry Smith PetscFunctionReturn(0); 68619717074SBarry Smith } 6878f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 68862cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 689a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6908f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 691a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6928f99978dSLois Curfman McInnes } 6938f99978dSLois Curfman McInnes } 694d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6953a40ed3dSBarry Smith PetscFunctionReturn(0); 6965e76c431SBarry Smith } 69704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6984a2ae208SSatish Balay #undef __FUNCT__ 69917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 7004b828684SBarry Smith /*@C 70117bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 7025e76c431SBarry Smith 703c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 704c7afd0dbSLois Curfman McInnes 7055e42470aSBarry Smith Input Parameters: 706c7afd0dbSLois Curfman McInnes + snes - the SNES context 707af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 7085e42470aSBarry Smith . x - current iterate 7095e42470aSBarry Smith . f - residual evaluated at x 7103c632250SBarry Smith . y - search direction 7115e42470aSBarry Smith . w - work vector 712dc357ad2SBarry Smith . fnorm - 2-norm of f 713dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 7145e42470aSBarry Smith 715c4a48953SLois Curfman McInnes Output Parameters: 7167f3332b4SBarry Smith + g - residual evaluated at new iterate w 717e106eecfSBarry Smith . w - new iterate (x + lambda*y) 7185e42470aSBarry Smith . gnorm - 2-norm of g 7195e42470aSBarry Smith . ynorm - 2-norm of search length 720a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 721fee21e36SBarry Smith 722ce9499c7SBarry Smith Options Database Keys: 723ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 724ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 725e106eecfSBarry 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) 726e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 727e106eecfSBarry Smith 7285e42470aSBarry Smith Notes: 72917bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7305e42470aSBarry Smith 73136851e7fSLois Curfman McInnes Level: advanced 73236851e7fSLois Curfman McInnes 733f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 734f59ffdedSLois Curfman McInnes 73517bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7365e42470aSBarry Smith @*/ 7377087cfbeSBarry Smith PetscErrorCode SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 7385e76c431SBarry Smith { 739406556e6SLois Curfman McInnes /* 740406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 741406556e6SLois Curfman McInnes minimization problem: 742406556e6SLois Curfman McInnes min z(x): R^n -> R, 743406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 744406556e6SLois Curfman McInnes */ 745e106eecfSBarry Smith PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 746aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 747f5b6597dSBarry Smith PetscScalar cinitslope; 7486b5873e3SBarry Smith #endif 749dfbe8321SBarry Smith PetscErrorCode ierr; 750a7cc72afSBarry Smith PetscInt count; 7514b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 752ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7535e76c431SBarry Smith 7543a40ed3dSBarry Smith PetscFunctionBegin; 755d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 756a7cc72afSBarry Smith *flag = PETSC_TRUE; 7575e76c431SBarry Smith 7583505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 75963b9fa5eSBarry Smith if (*ynorm == 0.0) { 76094298653SBarry Smith if (neP->monitor) { 76194298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 76294298653SBarry Smith } 763b37302c6SLois Curfman McInnes *gnorm = fnorm; 764e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 765b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 766a7cc72afSBarry Smith *flag = PETSC_FALSE; 767ad922ac9SBarry Smith goto theend2; 76894a9d846SBarry Smith } 769e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 770e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 771e106eecfSBarry Smith *ynorm = neP->maxstep; 7725e76c431SBarry Smith } 773ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 774e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 775a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 776aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 777a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 778329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7795e42470aSBarry Smith #else 780a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7815e42470aSBarry Smith #endif 7825e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7835e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 784f8833479SBarry Smith ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 7855e42470aSBarry Smith 786e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 78743e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 788ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 78943e71028SBarry Smith *flag = PETSC_FALSE; 79043e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 79143e71028SBarry Smith goto theend2; 79243e71028SBarry Smith } 7934936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 7944936397dSBarry Smith if (snes->domainerror) { 7954936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7964936397dSBarry Smith PetscFunctionReturn(0); 79719717074SBarry Smith } 798cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 79962cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 800e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 80194298653SBarry Smith if (neP->monitor) { 80294298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 80394298653SBarry Smith } 804ad922ac9SBarry Smith goto theend2; 8055e42470aSBarry Smith } 8065e42470aSBarry Smith 8075e42470aSBarry Smith /* Fit points with quadratic */ 808313b5538SBarry Smith lambda = 1.0; 8095e42470aSBarry Smith count = 1; 8105ca10a37SBarry Smith while (PETSC_TRUE) { 8115e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 81294298653SBarry Smith if (neP->monitor) { 81394298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 81494298653SBarry 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); 81594298653SBarry Smith } 816e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 81743e71028SBarry Smith *flag = PETSC_FALSE; 81843e71028SBarry Smith break; 8195e42470aSBarry Smith } 820a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 821329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 822329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 823329f5518SBarry Smith else lambda = lambdatemp; 824275d25c3SBarry Smith 825e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 82643e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 827ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 828d9822059SBarry Smith ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 829ed98166eSMatthew Knepley *flag = PETSC_FALSE; 83043e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 831ed98166eSMatthew Knepley break; 832ed98166eSMatthew Knepley } 833f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8344936397dSBarry Smith if (snes->domainerror) { 8354936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8364936397dSBarry Smith PetscFunctionReturn(0); 83719717074SBarry Smith } 838cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 83962cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 840e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 84194298653SBarry Smith if (neP->monitor) { 84294298653SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 84394298653SBarry Smith } 84406259719SBarry Smith break; 8455e42470aSBarry Smith } 8465e42470aSBarry Smith count++; 8475e42470aSBarry Smith } 848ad922ac9SBarry Smith theend2: 8498f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8503c632250SBarry Smith if (neP->postcheckstep) { 8513c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8523c632250SBarry Smith if (changed_y) { 85379f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8543c632250SBarry Smith } 8553c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8563c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 8574936397dSBarry Smith if (snes->domainerror) { 8584936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8594936397dSBarry Smith PetscFunctionReturn(0); 86019717074SBarry Smith } 8618f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 862a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8638f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 864a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 86562cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8668f99978dSLois Curfman McInnes } 8678f99978dSLois Curfman McInnes } 868d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8693a40ed3dSBarry Smith PetscFunctionReturn(0); 8705e76c431SBarry Smith } 8712343ba6eSBarry Smith 87204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8734a2ae208SSatish Balay #undef __FUNCT__ 87417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 875c9e6a524SLois Curfman McInnes /*@C 87617bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8774b27c08aSLois Curfman McInnes by the method SNESLS. 8785e76c431SBarry Smith 8795e76c431SBarry Smith Input Parameters: 880c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 881af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 882c7afd0dbSLois Curfman McInnes - func - pointer to int function 8835e76c431SBarry Smith 8843f9fe445SBarry Smith Logically Collective on SNES 885fee21e36SBarry Smith 886c4a48953SLois Curfman McInnes Available Routines: 88717bae607SBarry Smith + SNESLineSearchCubic() - default line search 88817bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 88917bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 89017bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8915e76c431SBarry Smith 892c4a48953SLois Curfman McInnes Options Database Keys: 8934b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8944b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 895e106eecfSBarry 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) 896e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 897c4a48953SLois Curfman McInnes 8985e76c431SBarry Smith Calling sequence of func: 899bd208895SLois Curfman McInnes .vb 900ace3abfcSBarry Smith func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 901bd208895SLois Curfman McInnes .ve 9025e76c431SBarry Smith 9035e76c431SBarry Smith Input parameters for func: 904c7afd0dbSLois Curfman McInnes + snes - nonlinear context 905af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 9065e76c431SBarry Smith . x - current iterate 9075e76c431SBarry Smith . f - residual evaluated at x 9083c632250SBarry Smith . y - search direction 909c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9105e76c431SBarry Smith 9115e76c431SBarry Smith Output parameters for func: 912c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9133c632250SBarry Smith . w - new iterate 9145e76c431SBarry Smith . gnorm - 2-norm of g 9155e76c431SBarry Smith . ynorm - 2-norm of search length 916a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 917f59ffdedSLois Curfman McInnes 91836851e7fSLois Curfman McInnes Level: advanced 91936851e7fSLois Curfman McInnes 920f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 921f59ffdedSLois Curfman McInnes 92217bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 92317bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9245e76c431SBarry Smith @*/ 9257087cfbeSBarry Smith PetscErrorCode SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx) 9265e76c431SBarry Smith { 9274ac538c5SBarry Smith PetscErrorCode ierr; 92882bf6240SBarry Smith 9293a40ed3dSBarry Smith PetscFunctionBegin; 9304ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSet_C",(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void*),(snes,func,lsctx));CHKERRQ(ierr); 9313a40ed3dSBarry Smith PetscFunctionReturn(0); 9325e76c431SBarry Smith } 9338e019c35SBarry Smith 934ace3abfcSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 93504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 936fb2e594dSBarry Smith EXTERN_C_BEGIN 9374a2ae208SSatish Balay #undef __FUNCT__ 93817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 9397087cfbeSBarry Smith PetscErrorCode SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 94082bf6240SBarry Smith { 94182bf6240SBarry Smith PetscFunctionBegin; 9424b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9434b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 94482bf6240SBarry Smith PetscFunctionReturn(0); 94582bf6240SBarry Smith } 946fb2e594dSBarry Smith EXTERN_C_END 94704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9484a2ae208SSatish Balay #undef __FUNCT__ 9493c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 950c8dd0c0dSLois Curfman McInnes /*@C 9513c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9524b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 953c8dd0c0dSLois Curfman McInnes 954c8dd0c0dSLois Curfman McInnes Input Parameters: 955c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9563c632250SBarry Smith . func - pointer to function 957c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 958c8dd0c0dSLois Curfman McInnes 9593f9fe445SBarry Smith Logically Collective on SNES 960c8dd0c0dSLois Curfman McInnes 961c8dd0c0dSLois Curfman McInnes Calling sequence of func: 962c8dd0c0dSLois Curfman McInnes .vb 963ace3abfcSBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool *changed_y,PetscBool *changed_w) 964c8dd0c0dSLois Curfman McInnes .ve 965b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 966b1ae27eaSLois Curfman McInnes on failure. 967c8dd0c0dSLois Curfman McInnes 968c8dd0c0dSLois Curfman McInnes Input parameters for func: 969c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 970c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9713c632250SBarry Smith . x - previous iterate 9723c632250SBarry Smith . y - new search direction and length 9733c632250SBarry Smith - w - current candidate iterate 974c8dd0c0dSLois Curfman McInnes 975c8dd0c0dSLois Curfman McInnes Output parameters for func: 9763c632250SBarry Smith + y - search direction (possibly changed) 9773c632250SBarry Smith . w - current iterate (possibly modified) 9783c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 9793c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 980c8dd0c0dSLois Curfman McInnes 981c8dd0c0dSLois Curfman McInnes Level: advanced 982c8dd0c0dSLois Curfman McInnes 9839e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 9849e247f21SBarry Smith 9853c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 9863c632250SBarry Smith 9873c632250SBarry Smith On input w = x + y 9883c632250SBarry Smith 98917bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 990b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 991ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9928f99978dSLois Curfman McInnes 99317bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 994ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 995ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 996ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 9979e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 998b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9998f99978dSLois Curfman McInnes 1000c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 1001c8dd0c0dSLois Curfman McInnes 100217bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1003c8dd0c0dSLois Curfman McInnes @*/ 10047087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx) 1005c8dd0c0dSLois Curfman McInnes { 10064ac538c5SBarry Smith PetscErrorCode ierr; 1007c8dd0c0dSLois Curfman McInnes 1008c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10094ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 1010c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1011c8dd0c0dSLois Curfman McInnes } 101294298653SBarry Smith 10139c3ca13aSBarry Smith #undef __FUNCT__ 10149c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10159c3ca13aSBarry Smith /*@C 10169c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10177e4bb74cSBarry Smith before the line search is called. 10189c3ca13aSBarry Smith 10199c3ca13aSBarry Smith Input Parameters: 10209c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10219c3ca13aSBarry Smith . func - pointer to function 10229c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10239c3ca13aSBarry Smith 10243f9fe445SBarry Smith Logically Collective on SNES 10259c3ca13aSBarry Smith 10269c3ca13aSBarry Smith Calling sequence of func: 10279c3ca13aSBarry Smith .vb 1028ace3abfcSBarry Smith int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool *changed_y) 10299c3ca13aSBarry Smith .ve 10309c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10319c3ca13aSBarry Smith on failure. 10329c3ca13aSBarry Smith 10339c3ca13aSBarry Smith Input parameters for func: 10349c3ca13aSBarry Smith + snes - nonlinear context 10359c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10369c3ca13aSBarry Smith . x - previous iterate 10379c3ca13aSBarry Smith - y - new search direction and length 10389c3ca13aSBarry Smith 10399c3ca13aSBarry Smith Output parameters for func: 10409c3ca13aSBarry Smith + y - search direction (possibly changed) 10419c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10429c3ca13aSBarry Smith 10439c3ca13aSBarry Smith Level: advanced 10449c3ca13aSBarry Smith 10459c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10469c3ca13aSBarry Smith 10479c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10489c3ca13aSBarry Smith 10497e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 10509c3ca13aSBarry Smith @*/ 10517087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx) 10529c3ca13aSBarry Smith { 10534ac538c5SBarry Smith PetscErrorCode ierr; 10549c3ca13aSBarry Smith 10559c3ca13aSBarry Smith PetscFunctionBegin; 10564ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 10579c3ca13aSBarry Smith PetscFunctionReturn(0); 10589c3ca13aSBarry Smith } 10599c3ca13aSBarry Smith 106094298653SBarry Smith #undef __FUNCT__ 106194298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor" 106294298653SBarry Smith /*@C 106394298653SBarry Smith SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search 106494298653SBarry Smith 106594298653SBarry Smith Input Parameters: 106694298653SBarry Smith + snes - nonlinear context obtained from SNESCreate() 106794298653SBarry Smith - flg - PETSC_TRUE to monitor the line search 106894298653SBarry Smith 10693f9fe445SBarry Smith Logically Collective on SNES 107094298653SBarry Smith 107194298653SBarry Smith Options Database: 107294298653SBarry Smith . -snes_ls_monitor 107394298653SBarry Smith 107494298653SBarry Smith Level: intermediate 107594298653SBarry Smith 107694298653SBarry Smith 107794298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 107894298653SBarry Smith @*/ 10797087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetMonitor(SNES snes,PetscBool flg) 108094298653SBarry Smith { 10814ac538c5SBarry Smith PetscErrorCode ierr; 108294298653SBarry Smith 108394298653SBarry Smith PetscFunctionBegin; 10844ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr); 108594298653SBarry Smith PetscFunctionReturn(0); 108694298653SBarry Smith } 108794298653SBarry Smith 1088c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1089ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/ 1090ace3abfcSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *); /* force argument to next function to not be extern C*/ 1091c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 10924a2ae208SSatish Balay #undef __FUNCT__ 10933c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 10947087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1095c8dd0c0dSLois Curfman McInnes { 1096c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10973c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 10983c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1099c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1100c8dd0c0dSLois Curfman McInnes } 1101c8dd0c0dSLois Curfman McInnes EXTERN_C_END 11029c3ca13aSBarry Smith 11039c3ca13aSBarry Smith EXTERN_C_BEGIN 11049c3ca13aSBarry Smith #undef __FUNCT__ 11059c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 11067087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 11079c3ca13aSBarry Smith { 11089c3ca13aSBarry Smith PetscFunctionBegin; 11099c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 11109c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 11119c3ca13aSBarry Smith PetscFunctionReturn(0); 11129c3ca13aSBarry Smith } 11139c3ca13aSBarry Smith EXTERN_C_END 1114329e7e40SMatthew Knepley 111594298653SBarry Smith EXTERN_C_BEGIN 111694298653SBarry Smith #undef __FUNCT__ 111794298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS" 11187087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetMonitor_LS(SNES snes,PetscBool flg) 111994298653SBarry Smith { 112094298653SBarry Smith SNES_LS *ls = (SNES_LS*)snes->data; 112194298653SBarry Smith PetscErrorCode ierr; 112294298653SBarry Smith 112394298653SBarry Smith PetscFunctionBegin; 112494298653SBarry Smith if (flg && !ls->monitor) { 112594298653SBarry Smith ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&ls->monitor);CHKERRQ(ierr); 112694298653SBarry Smith } else if (!flg && ls->monitor) { 11276bf464f9SBarry Smith ierr = PetscViewerASCIIMonitorDestroy(&ls->monitor);CHKERRQ(ierr); 112894298653SBarry Smith } 112994298653SBarry Smith PetscFunctionReturn(0); 113094298653SBarry Smith } 113194298653SBarry Smith EXTERN_C_END 113294298653SBarry Smith 1133329e7e40SMatthew Knepley /* 11344b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 113504d965bbSLois Curfman McInnes 113604d965bbSLois Curfman McInnes Input Parameters: 113704d965bbSLois Curfman McInnes . SNES - the SNES context 113804d965bbSLois Curfman McInnes . viewer - visualization context 113904d965bbSLois Curfman McInnes 114004d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 114104d965bbSLois Curfman McInnes */ 11424a2ae208SSatish Balay #undef __FUNCT__ 11434b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11446849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1145a935fc98SLois Curfman McInnes { 11464b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11472fc52814SBarry Smith const char *cstr; 1148dfbe8321SBarry Smith PetscErrorCode ierr; 1149ace3abfcSBarry Smith PetscBool iascii; 1150a935fc98SLois Curfman McInnes 11513a40ed3dSBarry Smith PetscFunctionBegin; 11522692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 115332077d6dSBarry Smith if (iascii) { 115417bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 115517bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 115617bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 115719bcc07fSBarry Smith else cstr = "unknown"; 1158b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1159e106eecfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 11605cd90555SBarry Smith } else { 1161e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 116219bcc07fSBarry Smith } 11633a40ed3dSBarry Smith PetscFunctionReturn(0); 1164a935fc98SLois Curfman McInnes } 116504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 116604d965bbSLois Curfman McInnes /* 11674b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 116804d965bbSLois Curfman McInnes 116904d965bbSLois Curfman McInnes Input Parameter: 117004d965bbSLois Curfman McInnes . snes - the SNES context 117104d965bbSLois Curfman McInnes 117204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 117304d965bbSLois Curfman McInnes */ 11744a2ae208SSatish Balay #undef __FUNCT__ 11754b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11766849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11775e42470aSBarry Smith { 11784b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1179e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1180dfbe8321SBarry Smith PetscErrorCode ierr; 1181a7cc72afSBarry Smith PetscInt indx; 1182ace3abfcSBarry Smith PetscBool flg,set; 11835e42470aSBarry Smith 11843a40ed3dSBarry Smith PetscFunctionBegin; 1185b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 11864b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 11874b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1188e106eecfSBarry Smith ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1189acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 119094298653SBarry Smith if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1191186905e3SBarry Smith 119217bae607SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 119325cdf11fSBarry Smith if (flg) { 119422e36657SBarry Smith switch (indx) { 1195b49fd9e1SBarry Smith case 0: 119617bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1197b49fd9e1SBarry Smith break; 1198b49fd9e1SBarry Smith case 1: 119917bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1200b49fd9e1SBarry Smith break; 1201b49fd9e1SBarry Smith case 2: 120217bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1203b49fd9e1SBarry Smith break; 1204b49fd9e1SBarry Smith case 3: 120517bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1206b49fd9e1SBarry Smith break; 12075e42470aSBarry Smith } 12085e42470aSBarry Smith } 1209b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 12103a40ed3dSBarry Smith PetscFunctionReturn(0); 12115e42470aSBarry Smith } 12124827ddcaSBarry Smith 12134827ddcaSBarry Smith EXTERN_C_BEGIN 12144827ddcaSBarry Smith extern PetscErrorCode SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal); 12154827ddcaSBarry Smith EXTERN_C_END 12164827ddcaSBarry Smith 121704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1218ebe3b25bSBarry Smith /*MC 1219ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 122004d965bbSLois Curfman McInnes 1221ebe3b25bSBarry Smith Options Database: 1222ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1223ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1224e106eecfSBarry 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) 1225acbee50cSBarry Smith . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1226acbee50cSBarry Smith - -snes_ls_monitor - print information about progress of line searches 1227acbee50cSBarry Smith 122804d965bbSLois Curfman McInnes 1229ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 123004d965bbSLois Curfman McInnes 1231ee3001cbSBarry Smith Level: beginner 1232ee3001cbSBarry Smith 123317bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 123417bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1235b3dd4ab5SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1236ebe3b25bSBarry Smith 1237ebe3b25bSBarry Smith M*/ 1238fb2e594dSBarry Smith EXTERN_C_BEGIN 12394a2ae208SSatish Balay #undef __FUNCT__ 12404b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 12417087cfbeSBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 12425e42470aSBarry Smith { 1243dfbe8321SBarry Smith PetscErrorCode ierr; 12444b27c08aSLois Curfman McInnes SNES_LS *neP; 12455e42470aSBarry Smith 12463a40ed3dSBarry Smith PetscFunctionBegin; 1247e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 1248e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 1249e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 1250e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 1251e7788613SBarry Smith snes->ops->view = SNESView_LS; 12526b8b9a38SLisandro Dalcin snes->ops->reset = SNESReset_LS; 12535e42470aSBarry Smith 125438f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 12555e42470aSBarry Smith snes->data = (void*)neP; 12565e42470aSBarry Smith neP->alpha = 1.e-4; 12575e42470aSBarry Smith neP->maxstep = 1.e8; 1258e106eecfSBarry Smith neP->minlambda = 1.e-12; 125917bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1260c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12613c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12623c632250SBarry Smith neP->postcheck = PETSC_NULL; 12633c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12643c632250SBarry Smith neP->precheck = PETSC_NULL; 126582bf6240SBarry Smith 126694298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr); 12674827ddcaSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr); 126894298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 126994298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 127094298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 127182bf6240SBarry Smith 12723a40ed3dSBarry Smith PetscFunctionReturn(0); 12735e42470aSBarry Smith } 1274fb2e594dSBarry Smith EXTERN_C_END 127582bf6240SBarry Smith 127682bf6240SBarry Smith 127782bf6240SBarry Smith 1278