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); 258f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(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); 398f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)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) { 678f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(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 { 1336849ba73SBarry Smith PetscErrorCode ierr; 134a7cc72afSBarry Smith PetscInt maxits,i,lits; 135ace3abfcSBarry Smith PetscBool lssucceed; 136112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 13785385478SLisandro Dalcin PetscReal fnorm,gnorm,xnorm=0,ynorm; 13885385478SLisandro Dalcin Vec Y,X,F,G,W; 1393d4c4710SBarry Smith KSPConvergedReason kspreason; 1405e76c431SBarry Smith 1413a40ed3dSBarry Smith PetscFunctionBegin; 1423d4c4710SBarry Smith snes->numFailures = 0; 1433d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 144184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 145184914b5SBarry Smith 1465e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1475e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1495e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1505e42470aSBarry Smith G = snes->work[1]; 1515e42470aSBarry Smith W = snes->work[2]; 1525e76c431SBarry Smith 1534c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1544c49b128SBarry Smith snes->iter = 0; 15575567043SBarry Smith snes->norm = 0.0; 1564c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15719717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1584936397dSBarry Smith if (snes->domainerror) { 1594936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1604936397dSBarry Smith PetscFunctionReturn(0); 1614936397dSBarry Smith } 1622613ca53SBarry Smith ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1632613ca53SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1642613ca53SBarry Smith ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 1652613ca53SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 16662cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1670f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1685e42470aSBarry Smith snes->norm = fnorm; 1690f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 17084c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 1717a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1723f149594SLisandro Dalcin 1733f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1743f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 17585385478SLisandro Dalcin /* test convergence */ 17685385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 17706ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 178d034289bSBarry Smith 1795e76c431SBarry Smith for (i=0; i<maxits; i++) { 1805e76c431SBarry Smith 181329e7e40SMatthew Knepley /* Call general purpose update function */ 182e7788613SBarry Smith if (snes->ops->update) { 183e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 184de8cb200SMatthew Knepley } 185329e7e40SMatthew Knepley 186ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1875334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 18894b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 18971f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 1903d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1913d4c4710SBarry Smith if (kspreason < 0) { 1923d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 193ef998cc9SBarry 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); 1943d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 195ab81109fSSatish Balay break; 1963d4c4710SBarry Smith } 1973d4c4710SBarry Smith } 1983d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1993f149594SLisandro Dalcin snes->linear_its += lits; 2003f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 20174637425SBarry Smith 202*ea630c6eSPeter Brune if (snes->ops->precheckstep) { 203ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 204*ea630c6eSPeter Brune ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr); 2059c3ca13aSBarry Smith } 2069c3ca13aSBarry Smith 207b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 2081e2582c4SBarry Smith ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 20985471664SBarry Smith } 21074637425SBarry Smith 211ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 212ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 213e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 214ea4d3ed3SLois Curfman McInnes */ 21585385478SLisandro Dalcin ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 2163f149594SLisandro Dalcin ynorm = 1; gnorm = fnorm; 217*ea630c6eSPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 2188f1a2a5eSBarry Smith ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 2194a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2204936397dSBarry Smith if (snes->domainerror) { 2214936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2224936397dSBarry Smith PetscFunctionReturn(0); 2234936397dSBarry Smith } 224a7cc72afSBarry Smith if (!lssucceed) { 22550ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 226ace3abfcSBarry Smith PetscBool ismin; 227647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2281e2582c4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 2298a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2303505fcc1SBarry Smith break; 2313505fcc1SBarry Smith } 23250ffb88aSMatthew Knepley } 23385385478SLisandro Dalcin /* Update function and solution vectors */ 23485385478SLisandro Dalcin fnorm = gnorm; 23585385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 23685385478SLisandro Dalcin ierr = VecCopy(W,X);CHKERRQ(ierr); 23785385478SLisandro Dalcin /* Monitor convergence */ 23885385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 23985385478SLisandro Dalcin snes->iter = i+1; 24085385478SLisandro Dalcin snes->norm = fnorm; 24185385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24285385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 2437a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 24485385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 24585385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 246e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2473f149594SLisandro Dalcin if (snes->reason) break; 24816c95cacSBarry Smith } 24952392280SLois Curfman McInnes if (i == maxits) { 250ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25185385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25252392280SLois Curfman McInnes } 2533a40ed3dSBarry Smith PetscFunctionReturn(0); 2545e76c431SBarry Smith } 25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25604d965bbSLois Curfman McInnes /* 2574b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2584b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 25904d965bbSLois Curfman McInnes 26004d965bbSLois Curfman McInnes Input Parameter: 26104d965bbSLois Curfman McInnes . snes - the SNES context 26204d965bbSLois Curfman McInnes . x - the solution vector 26304d965bbSLois Curfman McInnes 26404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26504d965bbSLois Curfman McInnes 26604d965bbSLois Curfman McInnes Notes: 2674b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 26804d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 26904d965bbSLois Curfman McInnes the call to SNESSolve(). 27004d965bbSLois Curfman McInnes */ 2714a2ae208SSatish Balay #undef __FUNCT__ 2724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 273dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2745e76c431SBarry Smith { 275dfbe8321SBarry Smith PetscErrorCode ierr; 2763a40ed3dSBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 27858c9b817SLisandro Dalcin ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 2805e76c431SBarry Smith } 28104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2826b8b9a38SLisandro Dalcin 2836b8b9a38SLisandro Dalcin #undef __FUNCT__ 2846b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS" 2856b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes) 2866b8b9a38SLisandro Dalcin { 2876b8b9a38SLisandro Dalcin PetscErrorCode ierr; 2886b8b9a38SLisandro Dalcin 2896b8b9a38SLisandro Dalcin PetscFunctionBegin; 2906b8b9a38SLisandro Dalcin if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2916b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2926b8b9a38SLisandro Dalcin } 2936b8b9a38SLisandro Dalcin 29404d965bbSLois Curfman McInnes /* 2954b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2964b27c08aSLois Curfman McInnes with SNESCreate_LS(). 29704d965bbSLois Curfman McInnes 29804d965bbSLois Curfman McInnes Input Parameter: 29904d965bbSLois Curfman McInnes . snes - the SNES context 30004d965bbSLois Curfman McInnes 301de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30204d965bbSLois Curfman McInnes */ 3034a2ae208SSatish Balay #undef __FUNCT__ 3044b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 305dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3065e76c431SBarry Smith { 307dfbe8321SBarry Smith PetscErrorCode ierr; 3083a40ed3dSBarry Smith 3093a40ed3dSBarry Smith PetscFunctionBegin; 3106b8b9a38SLisandro Dalcin ierr = SNESReset_LS(snes);CHKERRQ(ierr); 3115c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 312557d3f75SLisandro Dalcin 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 3145e76c431SBarry Smith } 31504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3164a2ae208SSatish Balay #undef __FUNCT__ 317*ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchNo_LS" 31882bf6240SBarry Smith 3194b828684SBarry Smith /*@C 32017bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3215e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3225e76c431SBarry Smith to serve as a template and is not recommended for general use. 3235e76c431SBarry Smith 3243f9fe445SBarry Smith Logically Collective on SNES and Vec 325c7afd0dbSLois Curfman McInnes 3265e76c431SBarry Smith Input Parameters: 327c7afd0dbSLois Curfman McInnes + snes - nonlinear context 328af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3295e76c431SBarry Smith . x - current iterate 3305e76c431SBarry Smith . f - residual evaluated at x 3313c632250SBarry Smith . y - search direction 332dc357ad2SBarry Smith . fnorm - 2-norm of f 333dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 3345e76c431SBarry Smith 335c4a48953SLois Curfman McInnes Output Parameters: 336c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3373c632250SBarry Smith . w - new iterate 3385e76c431SBarry Smith . gnorm - 2-norm of g 3395e76c431SBarry Smith . ynorm - 2-norm of search length 340a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 341fee21e36SBarry Smith 342c4a48953SLois Curfman McInnes Options Database Key: 34317bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 344c4a48953SLois Curfman McInnes 34536851e7fSLois Curfman McInnes Level: advanced 34636851e7fSLois Curfman McInnes 34728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 34828ae5a21SLois Curfman McInnes 34917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 35017bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3515e76c431SBarry Smith @*/ 352*ea630c6eSPeter Brune PetscErrorCode SNESLineSearchNo_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 3535e76c431SBarry Smith { 354dfbe8321SBarry Smith PetscErrorCode ierr; 355ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3565334005bSBarry Smith 3573a40ed3dSBarry Smith PetscFunctionBegin; 358a7cc72afSBarry Smith *flag = PETSC_TRUE; 359d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 360a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 361*ea630c6eSPeter Brune ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr); /* w <- x - y */ 362*ea630c6eSPeter Brune if (snes->ops->postcheckstep) { 363*ea630c6eSPeter Brune ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3648f99978dSLois Curfman McInnes } 3653c632250SBarry Smith if (changed_y) { 366*ea630c6eSPeter Brune ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr); /* w <- x - y */ 3673c632250SBarry Smith } 3684936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 3694936397dSBarry Smith if (!snes->domainerror) { 370a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 37162cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3724936397dSBarry Smith } 373d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3743a40ed3dSBarry Smith PetscFunctionReturn(0); 3755e76c431SBarry Smith } 37604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3774a2ae208SSatish Balay #undef __FUNCT__ 378*ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchNoNorms_LS" 37982bf6240SBarry Smith 38029e0b56fSBarry Smith /*@C 38117bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 38229e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 38329e0b56fSBarry Smith even compute the norm of the function or search direction; this 38429e0b56fSBarry Smith is intended only when you know the full step is fine and are 38529e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 386c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 387c7afd0dbSLois Curfman McInnes 3883f9fe445SBarry Smith Logically Collective on SNES and Vec 38929e0b56fSBarry Smith 39029e0b56fSBarry Smith Input Parameters: 391c7afd0dbSLois Curfman McInnes + snes - nonlinear context 392af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 39329e0b56fSBarry Smith . x - current iterate 39429e0b56fSBarry Smith . f - residual evaluated at x 3953c632250SBarry Smith . y - search direction 396dc357ad2SBarry Smith . fnorm - 2-norm of f 397dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 39829e0b56fSBarry Smith 39929e0b56fSBarry Smith Output Parameters: 400c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4013c632250SBarry Smith . w - new iterate 40229e0b56fSBarry Smith . gnorm - not changed 40329e0b56fSBarry Smith . ynorm - not changed 404a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 405fee21e36SBarry Smith 40629e0b56fSBarry Smith Options Database Key: 40717bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 40829e0b56fSBarry Smith 4098cbba510SBarry Smith Notes: 41017bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 411ea56f5baSLois Curfman McInnes either the options 412ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 413ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 414329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 415329f5518SBarry Smith otherwise, the SNES solver will generate an error. 416329f5518SBarry Smith 417329f5518SBarry Smith During the final iteration this will not evaluate the function at 418329f5518SBarry Smith the solution point. This is to save a function evaluation while 419329f5518SBarry Smith using pseudo-timestepping. 4208cbba510SBarry Smith 421ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 422a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 423ea56f5baSLois Curfman McInnes correct, since they are not computed. 424ea56f5baSLois Curfman McInnes 425ea56f5baSLois Curfman McInnes Level: advanced 4268cbba510SBarry Smith 42729e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 42829e0b56fSBarry Smith 42917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 43017bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 43129e0b56fSBarry Smith @*/ 432*ea630c6eSPeter Brune PetscErrorCode SNESLineSearchNoNorms_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 43329e0b56fSBarry Smith { 434dfbe8321SBarry Smith PetscErrorCode ierr; 435ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 43629e0b56fSBarry Smith 4373a40ed3dSBarry Smith PetscFunctionBegin; 438a7cc72afSBarry Smith *flag = PETSC_TRUE; 439d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 440*ea630c6eSPeter Brune ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr); /* w <- x - y */ 441*ea630c6eSPeter Brune if (snes->ops->postcheckstep) { 442*ea630c6eSPeter Brune ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4433c632250SBarry Smith } 4443c632250SBarry Smith if (changed_y) { 445*ea630c6eSPeter Brune ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr); /* w <- x - y */ 4468f99978dSLois Curfman McInnes } 447329f5518SBarry Smith 448329f5518SBarry Smith /* don't evaluate function the last time through */ 449329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4504936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 451329f5518SBarry Smith } 452d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4533a40ed3dSBarry Smith PetscFunctionReturn(0); 45429e0b56fSBarry Smith } 45504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4564a2ae208SSatish Balay #undef __FUNCT__ 457*ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchCubic_LS" 4584b828684SBarry Smith /*@C 45917bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4605e76c431SBarry Smith 461c7afd0dbSLois Curfman McInnes Collective on SNES 462c7afd0dbSLois Curfman McInnes 4635e76c431SBarry Smith Input Parameters: 464c7afd0dbSLois Curfman McInnes + snes - nonlinear context 465af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4665e76c431SBarry Smith . x - current iterate 4675e76c431SBarry Smith . f - residual evaluated at x 4683c632250SBarry Smith . y - search direction 469dc357ad2SBarry Smith . fnorm - 2-norm of f 470dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 4715e76c431SBarry Smith 472393d2d9aSLois Curfman McInnes Output Parameters: 473c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4743c632250SBarry Smith . w - new iterate 4755e76c431SBarry Smith . gnorm - 2-norm of g 4765e76c431SBarry Smith . ynorm - 2-norm of search length 477a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 478fee21e36SBarry Smith 479c4a48953SLois Curfman McInnes Options Database Key: 480e106eecfSBarry Smith + -snes_ls cubic - Activates SNESLineSearchCubic() 481e106eecfSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 482e106eecfSBarry 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) 483e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 484e106eecfSBarry Smith 485c4a48953SLois Curfman McInnes 4865e76c431SBarry Smith Notes: 4875e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4885e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4895e76c431SBarry Smith 49036851e7fSLois Curfman McInnes Level: advanced 49136851e7fSLois Curfman McInnes 49228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 49328ae5a21SLois Curfman McInnes 49417bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 4955e76c431SBarry Smith @*/ 496*ea630c6eSPeter Brune PetscErrorCode SNESLineSearchCubic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 4975e76c431SBarry Smith { 498406556e6SLois Curfman McInnes /* 499406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 500406556e6SLois Curfman McInnes minimization problem: 501406556e6SLois Curfman McInnes min z(x): R^n -> R, 502406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 503406556e6SLois Curfman McInnes */ 504406556e6SLois Curfman McInnes 505e106eecfSBarry Smith PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 506e106eecfSBarry Smith PetscReal minlambda,lambda,lambdatemp; 507aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 508f5b6597dSBarry Smith PetscScalar cinitslope; 5096b5873e3SBarry Smith #endif 510dfbe8321SBarry Smith PetscErrorCode ierr; 511a7cc72afSBarry Smith PetscInt count; 512ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 51394298653SBarry Smith MPI_Comm comm; 5145e76c431SBarry Smith 5153a40ed3dSBarry Smith PetscFunctionBegin; 51694298653SBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 517d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 518a7cc72afSBarry Smith *flag = PETSC_TRUE; 5195e76c431SBarry Smith 520cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 52175567043SBarry Smith if (*ynorm == 0.0) { 522*ea630c6eSPeter Brune if (snes->ls_monitor) { 523*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 524*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 525*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 52694298653SBarry Smith } 527a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5283c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 529a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 530a7cc72afSBarry Smith *flag = PETSC_FALSE; 531ad922ac9SBarry Smith goto theend1; 53294a9d846SBarry Smith } 533*ea630c6eSPeter Brune if (*ynorm > snes->maxstep) { /* Step too big, so scale back */ 534*ea630c6eSPeter Brune if (snes->ls_monitor) { 535*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 536*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr); 537*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 53894298653SBarry Smith } 539*ea630c6eSPeter Brune ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr); 540*ea630c6eSPeter Brune *ynorm = snes->maxstep; 5415e76c431SBarry Smith } 542ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 543*ea630c6eSPeter Brune minlambda = snes->steptol/rellength; 544a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 545aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 546a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 547329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5485e42470aSBarry Smith #else 549a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5505e42470aSBarry Smith #endif 5515e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5525e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5535e76c431SBarry Smith 554e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 55543e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 556ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 55743e71028SBarry Smith *flag = PETSC_FALSE; 55843e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 55943e71028SBarry Smith goto theend1; 56043e71028SBarry Smith } 5614936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5624936397dSBarry Smith if (snes->domainerror) { 5634936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5644936397dSBarry Smith PetscFunctionReturn(0); 56519717074SBarry Smith } 566313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 56762cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5688f1a2a5eSBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr); 569*ea630c6eSPeter Brune if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */ 570*ea630c6eSPeter Brune if (snes->ls_monitor) { 571*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 572*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr); 573*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 57494298653SBarry Smith } 575ad922ac9SBarry Smith goto theend1; 5765e76c431SBarry Smith } 5775e76c431SBarry Smith 5785e76c431SBarry Smith /* Fit points with quadratic */ 579313b5538SBarry Smith lambda = 1.0; 580a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5815e76c431SBarry Smith lambdaprev = lambda; 5825e76c431SBarry Smith gnormprev = *gnorm; 583329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 584ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 585ddd12767SBarry Smith else lambda = lambdatemp; 586275d25c3SBarry Smith 587e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 58843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 589ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 59043e71028SBarry Smith *flag = PETSC_FALSE; 59143e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 59243e71028SBarry Smith goto theend1; 59343e71028SBarry Smith } 594f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5954936397dSBarry Smith if (snes->domainerror) { 5964936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5974936397dSBarry Smith PetscFunctionReturn(0); 59819717074SBarry Smith } 599cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 60062cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 601*ea630c6eSPeter Brune if (snes->ls_monitor) { 602*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 603*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr); 604*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 60594298653SBarry Smith } 606*ea630c6eSPeter Brune if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */ 607*ea630c6eSPeter Brune if (snes->ls_monitor) { 608*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 609*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 610*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 61194298653SBarry Smith } 612ad922ac9SBarry Smith goto theend1; 6135e76c431SBarry Smith } 6145e76c431SBarry Smith 6155e76c431SBarry Smith /* Fit points with cubic */ 6165e76c431SBarry Smith count = 1; 617bb9195b6SBarry Smith while (PETSC_TRUE) { 618dc357ad2SBarry Smith if (lambda <= minlambda) { 619*ea630c6eSPeter Brune if (snes->ls_monitor) { 620*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 621*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 622*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_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); 623*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 62494298653SBarry Smith } 62543e71028SBarry Smith *flag = PETSC_FALSE; 62643e71028SBarry Smith break; 6275e76c431SBarry Smith } 6285d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6295d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 630ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6312b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6325e76c431SBarry Smith d = b*b - 3*a*initslope; 6335e76c431SBarry Smith if (d < 0.0) d = 0.0; 6345e76c431SBarry Smith if (a == 0.0) { 6355e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6365e76c431SBarry Smith } else { 6375e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6385e76c431SBarry Smith } 6395e76c431SBarry Smith lambdaprev = lambda; 6405e76c431SBarry Smith gnormprev = *gnorm; 641329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 642329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6435e76c431SBarry Smith else lambda = lambdatemp; 644e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 64543e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 646ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 647ae15b995SBarry 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); 648ed98166eSMatthew Knepley *flag = PETSC_FALSE; 64943e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 650ed98166eSMatthew Knepley break; 651ed98166eSMatthew Knepley } 652f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6534936397dSBarry Smith if (snes->domainerror) { 6544936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6554936397dSBarry Smith PetscFunctionReturn(0); 65619717074SBarry Smith } 657cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 65862cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 659*ea630c6eSPeter Brune if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* is reduction enough? */ 660*ea630c6eSPeter Brune if (snes->ls_monitor) { 6618f1a2a5eSBarry Smith ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr); 66294298653SBarry Smith } 66343e71028SBarry Smith break; 6642b022350SLois Curfman McInnes } else { 665*ea630c6eSPeter Brune if (snes->ls_monitor) { 6668f1a2a5eSBarry Smith ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr); 66794298653SBarry Smith } 6685e76c431SBarry Smith } 6695e76c431SBarry Smith count++; 6705e76c431SBarry Smith } 671ad922ac9SBarry Smith theend1: 6728f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 673*ea630c6eSPeter Brune if (snes->ops->postcheckstep && *flag) { 674*ea630c6eSPeter Brune ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6753c632250SBarry Smith if (changed_y) { 67679f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6773c632250SBarry Smith } 6783c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 679f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6804936397dSBarry Smith if (snes->domainerror) { 6814936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6824936397dSBarry Smith PetscFunctionReturn(0); 68319717074SBarry Smith } 6848f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 68562cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 686a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6878f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 688a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6898f99978dSLois Curfman McInnes } 6908f99978dSLois Curfman McInnes } 691d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6923a40ed3dSBarry Smith PetscFunctionReturn(0); 6935e76c431SBarry Smith } 69404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6954a2ae208SSatish Balay #undef __FUNCT__ 696*ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_LS" 6974b828684SBarry Smith /*@C 698*ea630c6eSPeter Brune SNESLineSearchQuadratic_LS - Performs a quadratic line search. 6995e76c431SBarry Smith 700c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 701c7afd0dbSLois Curfman McInnes 7025e42470aSBarry Smith Input Parameters: 703c7afd0dbSLois Curfman McInnes + snes - the SNES context 704af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 7055e42470aSBarry Smith . x - current iterate 7065e42470aSBarry Smith . f - residual evaluated at x 7073c632250SBarry Smith . y - search direction 708dc357ad2SBarry Smith . fnorm - 2-norm of f 709dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 7105e42470aSBarry Smith 711c4a48953SLois Curfman McInnes Output Parameters: 7127f3332b4SBarry Smith + g - residual evaluated at new iterate w 713e106eecfSBarry Smith . w - new iterate (x + lambda*y) 7145e42470aSBarry Smith . gnorm - 2-norm of g 7155e42470aSBarry Smith . ynorm - 2-norm of search length 716a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 717fee21e36SBarry Smith 718ce9499c7SBarry Smith Options Database Keys: 719ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 720ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 721e106eecfSBarry 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) 722e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 723e106eecfSBarry Smith 7245e42470aSBarry Smith Notes: 72517bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7265e42470aSBarry Smith 72736851e7fSLois Curfman McInnes Level: advanced 72836851e7fSLois Curfman McInnes 729f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 730f59ffdedSLois Curfman McInnes 73117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7325e42470aSBarry Smith @*/ 733*ea630c6eSPeter Brune PetscErrorCode SNESLineSearchQuadratic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 7345e76c431SBarry Smith { 735406556e6SLois Curfman McInnes /* 736406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 737406556e6SLois Curfman McInnes minimization problem: 738406556e6SLois Curfman McInnes min z(x): R^n -> R, 739406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 740406556e6SLois Curfman McInnes */ 741e106eecfSBarry Smith PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 742aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 743f5b6597dSBarry Smith PetscScalar cinitslope; 7446b5873e3SBarry Smith #endif 745dfbe8321SBarry Smith PetscErrorCode ierr; 746a7cc72afSBarry Smith PetscInt count; 747ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7485e76c431SBarry Smith 7493a40ed3dSBarry Smith PetscFunctionBegin; 750d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 751a7cc72afSBarry Smith *flag = PETSC_TRUE; 7525e76c431SBarry Smith 7533505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 75463b9fa5eSBarry Smith if (*ynorm == 0.0) { 755*ea630c6eSPeter Brune if (snes->ls_monitor) { 756*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 757*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 758*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 75994298653SBarry Smith } 760b37302c6SLois Curfman McInnes *gnorm = fnorm; 761e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 762b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 763a7cc72afSBarry Smith *flag = PETSC_FALSE; 764ad922ac9SBarry Smith goto theend2; 76594a9d846SBarry Smith } 766*ea630c6eSPeter Brune if (*ynorm > snes->maxstep) { /* Step too big, so scale back */ 767*ea630c6eSPeter Brune ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr); 768*ea630c6eSPeter Brune *ynorm = snes->maxstep; 7695e76c431SBarry Smith } 770ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 771*ea630c6eSPeter Brune minlambda = snes->steptol/rellength; 772a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 773aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 774a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 775329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7765e42470aSBarry Smith #else 777a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7785e42470aSBarry Smith #endif 7795e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7805e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7818f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr); 7825e42470aSBarry Smith 783e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 78443e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 785ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 78643e71028SBarry Smith *flag = PETSC_FALSE; 78743e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 78843e71028SBarry Smith goto theend2; 78943e71028SBarry Smith } 7904936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 7914936397dSBarry Smith if (snes->domainerror) { 7924936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7934936397dSBarry Smith PetscFunctionReturn(0); 79419717074SBarry Smith } 795cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 79662cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 797*ea630c6eSPeter Brune if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */ 798*ea630c6eSPeter Brune if (snes->ls_monitor) { 799*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 800*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line Search: Using full step\n");CHKERRQ(ierr); 801*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 80294298653SBarry Smith } 803ad922ac9SBarry Smith goto theend2; 8045e42470aSBarry Smith } 8055e42470aSBarry Smith 8065e42470aSBarry Smith /* Fit points with quadratic */ 807313b5538SBarry Smith lambda = 1.0; 8085e42470aSBarry Smith count = 1; 8095ca10a37SBarry Smith while (PETSC_TRUE) { 8105e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 811*ea630c6eSPeter Brune if (snes->ls_monitor) { 812*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 813*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 814*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 815*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 81694298653SBarry Smith } 817e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 81843e71028SBarry Smith *flag = PETSC_FALSE; 81943e71028SBarry Smith break; 8205e42470aSBarry Smith } 821a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 822329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 823329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 824329f5518SBarry Smith else lambda = lambdatemp; 825275d25c3SBarry Smith 826e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 82743e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 828ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 829d9822059SBarry 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); 830ed98166eSMatthew Knepley *flag = PETSC_FALSE; 83143e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 832ed98166eSMatthew Knepley break; 833ed98166eSMatthew Knepley } 834f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8354936397dSBarry Smith if (snes->domainerror) { 8364936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8374936397dSBarry Smith PetscFunctionReturn(0); 83819717074SBarry Smith } 839cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 84062cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 841*ea630c6eSPeter Brune if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */ 842*ea630c6eSPeter Brune if (snes->ls_monitor) { 843*ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 844*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr); 845*ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 84694298653SBarry Smith } 84706259719SBarry Smith break; 8485e42470aSBarry Smith } 8495e42470aSBarry Smith count++; 8505e42470aSBarry Smith } 851ad922ac9SBarry Smith theend2: 8528f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 853*ea630c6eSPeter Brune if (snes->ops->postcheckstep) { 854*ea630c6eSPeter Brune ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8553c632250SBarry Smith if (changed_y) { 85679f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8573c632250SBarry Smith } 8583c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8593c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 8604936397dSBarry Smith if (snes->domainerror) { 8614936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8624936397dSBarry Smith PetscFunctionReturn(0); 86319717074SBarry Smith } 8648f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 865a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8668f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 867a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 86862cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8698f99978dSLois Curfman McInnes } 8708f99978dSLois Curfman McInnes } 871d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8723a40ed3dSBarry Smith PetscFunctionReturn(0); 8735e76c431SBarry Smith } 87404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 87594298653SBarry Smith 876329e7e40SMatthew Knepley /* 8774b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 87804d965bbSLois Curfman McInnes 87904d965bbSLois Curfman McInnes Input Parameters: 88004d965bbSLois Curfman McInnes . SNES - the SNES context 88104d965bbSLois Curfman McInnes . viewer - visualization context 88204d965bbSLois Curfman McInnes 88304d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 88404d965bbSLois Curfman McInnes */ 8854a2ae208SSatish Balay #undef __FUNCT__ 8864b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 8876849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 888a935fc98SLois Curfman McInnes { 8892fc52814SBarry Smith const char *cstr; 890dfbe8321SBarry Smith PetscErrorCode ierr; 891ace3abfcSBarry Smith PetscBool iascii; 892a935fc98SLois Curfman McInnes 8933a40ed3dSBarry Smith PetscFunctionBegin; 8942692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 89532077d6dSBarry Smith if (iascii) { 896*ea630c6eSPeter Brune if (snes->ops->linesearch == SNESLineSearchNo_LS) cstr = "SNESLineSearchNo"; 897*ea630c6eSPeter Brune if (snes->ops->linesearch == SNESLineSearchNoNorms_LS) cstr = "SNESLineSearchNoNorms"; 898*ea630c6eSPeter Brune else if (snes->ops->linesearch == SNESLineSearchQuadratic_LS) cstr = "SNESLineSearchQuadratic"; 899*ea630c6eSPeter Brune else if (snes->ops->linesearch == SNESLineSearchCubic_LS) cstr = "SNESLineSearchCubic"; 90019bcc07fSBarry Smith else cstr = "unknown"; 901b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 902*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)snes->ls_alpha,(double)snes->maxstep,(double)snes->steptol);CHKERRQ(ierr); 903*ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr); 90419bcc07fSBarry Smith } 9053a40ed3dSBarry Smith PetscFunctionReturn(0); 906a935fc98SLois Curfman McInnes } 90704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 90804d965bbSLois Curfman McInnes /* 9094b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 91004d965bbSLois Curfman McInnes 91104d965bbSLois Curfman McInnes Input Parameter: 91204d965bbSLois Curfman McInnes . snes - the SNES context 91304d965bbSLois Curfman McInnes 91404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 91504d965bbSLois Curfman McInnes */ 9164a2ae208SSatish Balay #undef __FUNCT__ 9174b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 9186849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 9195e42470aSBarry Smith { 920dfbe8321SBarry Smith PetscErrorCode ierr; 9215e42470aSBarry Smith 9223a40ed3dSBarry Smith PetscFunctionBegin; 923b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 924b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 9253a40ed3dSBarry Smith PetscFunctionReturn(0); 9265e42470aSBarry Smith } 9274827ddcaSBarry Smith 9284827ddcaSBarry Smith EXTERN_C_BEGIN 9294827ddcaSBarry Smith extern PetscErrorCode SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal); 9304827ddcaSBarry Smith EXTERN_C_END 9314827ddcaSBarry Smith 93204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 933ebe3b25bSBarry Smith /*MC 934ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 93504d965bbSLois Curfman McInnes 936ebe3b25bSBarry Smith Options Database: 937ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 93855d4ea13SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 939e106eecfSBarry 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) 940acbee50cSBarry Smith . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 94155d4ea13SBarry Smith . -snes_ls_monitor - print information about progress of line searches 94255d4ea13SBarry Smith - -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms 943acbee50cSBarry Smith 94404d965bbSLois Curfman McInnes 945ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 94604d965bbSLois Curfman McInnes 947ee3001cbSBarry Smith Level: beginner 948ee3001cbSBarry Smith 94917bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 95017bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 951b3dd4ab5SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 952ebe3b25bSBarry Smith 953ebe3b25bSBarry Smith M*/ 954fb2e594dSBarry Smith EXTERN_C_BEGIN 9554a2ae208SSatish Balay #undef __FUNCT__ 9564b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 9577087cfbeSBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 9585e42470aSBarry Smith { 959dfbe8321SBarry Smith PetscErrorCode ierr; 9604b27c08aSLois Curfman McInnes SNES_LS *neP; 9615e42470aSBarry Smith 9623a40ed3dSBarry Smith PetscFunctionBegin; 963e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 964e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 965e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 966e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 967e7788613SBarry Smith snes->ops->view = SNESView_LS; 9686b8b9a38SLisandro Dalcin snes->ops->reset = SNESReset_LS; 9695e42470aSBarry Smith 97042f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 97142f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 97238f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 9735e42470aSBarry Smith snes->data = (void*)neP; 974*ea630c6eSPeter Brune snes->ops->linesearch = SNESLineSearchCubic_LS; 975*ea630c6eSPeter Brune snes->ops->linesearchno = SNESLineSearchNo_LS; 976*ea630c6eSPeter Brune snes->ops->linesearchquadratic = SNESLineSearchQuadratic_LS; 977*ea630c6eSPeter Brune snes->ops->linesearchcubic = SNESLineSearchCubic_LS; 978*ea630c6eSPeter Brune snes->ops->linesearchexact = PETSC_NULL; 979*ea630c6eSPeter Brune snes->ops->linesearchtest = PETSC_NULL; 980*ea630c6eSPeter Brune snes->lsP = PETSC_NULL; 981*ea630c6eSPeter Brune snes->ops->postcheckstep = PETSC_NULL; 982*ea630c6eSPeter Brune snes->postcheck = PETSC_NULL; 983*ea630c6eSPeter Brune snes->ops->precheckstep = PETSC_NULL; 984*ea630c6eSPeter Brune snes->precheck = PETSC_NULL; 98582bf6240SBarry Smith 9863a40ed3dSBarry Smith PetscFunctionReturn(0); 9875e42470aSBarry Smith } 988fb2e594dSBarry Smith EXTERN_C_END 989