15e76c431SBarry Smith 2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h> 35e42470aSBarry Smith 4*490ca623SBarry Smith const char *SNESLineSearchTypes[] = {"BASIC","BASICNONORMS","QUADRATIC","CUBIC","SNESLineSearchType","SNES_LS_",0}; 5*490ca623SBarry Smith 68a5d9ee4SBarry Smith /* 720f52e01SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 820f52e01SBarry 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 936669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1020f52e01SBarry Smith for this trick. One assumes that the probability that W is in the null space of J is very, very small. 118a5d9ee4SBarry Smith */ 124a2ae208SSatish Balay #undef __FUNCT__ 134a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 14ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 158a5d9ee4SBarry Smith { 168a5d9ee4SBarry Smith PetscReal a1; 17dfbe8321SBarry Smith PetscErrorCode ierr; 18ace3abfcSBarry Smith PetscBool hastranspose; 198a5d9ee4SBarry Smith 208a5d9ee4SBarry Smith PetscFunctionBegin; 218a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2236669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2336669109SBarry Smith if (hastranspose) { 248a5d9ee4SBarry Smith /* Compute || J^T F|| */ 258a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 268a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 278f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 288a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2936669109SBarry Smith } else { 3036669109SBarry Smith Vec work; 31ea709b57SSatish Balay PetscScalar result; 3236669109SBarry Smith PetscReal wnorm; 3336669109SBarry Smith 34abb0e124SMatthew Knepley ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3636669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 396bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 4036669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 418f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 4236669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4336669109SBarry Smith } 448a5d9ee4SBarry Smith PetscFunctionReturn(0); 458a5d9ee4SBarry Smith } 468a5d9ee4SBarry Smith 4774637425SBarry Smith /* 485ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4974637425SBarry Smith */ 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 521e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 5374637425SBarry Smith { 5474637425SBarry Smith PetscReal a1,a2; 55dfbe8321SBarry Smith PetscErrorCode ierr; 56ace3abfcSBarry Smith PetscBool hastranspose; 5774637425SBarry Smith 5874637425SBarry Smith PetscFunctionBegin; 5974637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 6074637425SBarry Smith if (hastranspose) { 6174637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6279f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6374637425SBarry Smith 6474637425SBarry Smith /* Compute || J^T W|| */ 6574637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6774637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 6875567043SBarry Smith if (a1 != 0.0) { 698f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 7085471664SBarry Smith } 7174637425SBarry Smith } 7274637425SBarry Smith PetscFunctionReturn(0); 7374637425SBarry Smith } 7474637425SBarry Smith 7504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7604d965bbSLois Curfman McInnes 7704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7894b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8004d965bbSLois Curfman McInnes respectively. 8104d965bbSLois Curfman McInnes 82fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 85fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 884b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 894b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9304d965bbSLois Curfman McInnes for all nonlinear solvers. 9404d965bbSLois Curfman McInnes 9504d965bbSLois Curfman McInnes Another key routine is: 9604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9904d965bbSLois Curfman McInnes SNESSolve() if necessary. 10004d965bbSLois Curfman McInnes 10104d965bbSLois Curfman McInnes Additional basic routines are: 10204d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10304d965bbSLois Curfman McInnes have actually been used. 10404d965bbSLois Curfman McInnes These are called by application codes via the interface routines 105186905e3SBarry Smith SNESView(). 10604d965bbSLois Curfman McInnes 10704d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10804d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10904d965bbSLois Curfman McInnes above description applies to these categories also. 11004d965bbSLois Curfman McInnes 11104d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1125e42470aSBarry Smith /* 1134b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11404d965bbSLois Curfman McInnes method with a line search. 1155e76c431SBarry Smith 11604d965bbSLois Curfman McInnes Input Parameters: 11704d965bbSLois Curfman McInnes . snes - the SNES context 1185e76c431SBarry Smith 11904d965bbSLois Curfman McInnes Output Parameter: 120c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12104d965bbSLois Curfman McInnes 12204d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1235e76c431SBarry Smith 1245e76c431SBarry Smith Notes: 1255e76c431SBarry Smith This implements essentially a truncated Newton method with a 1265e76c431SBarry Smith line search. By default a cubic backtracking line search 1275e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1285e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 129393d2d9aSLois Curfman McInnes and Schnabel. 1305e76c431SBarry Smith */ 1314a2ae208SSatish Balay #undef __FUNCT__ 1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 133dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1345e76c431SBarry Smith { 1354b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1366849ba73SBarry Smith PetscErrorCode ierr; 137a7cc72afSBarry Smith PetscInt maxits,i,lits; 138ace3abfcSBarry Smith PetscBool lssucceed; 139112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 14085385478SLisandro Dalcin PetscReal fnorm,gnorm,xnorm=0,ynorm; 14185385478SLisandro Dalcin Vec Y,X,F,G,W; 1423d4c4710SBarry Smith KSPConvergedReason kspreason; 1435e76c431SBarry Smith 1443a40ed3dSBarry Smith PetscFunctionBegin; 1453d4c4710SBarry Smith snes->numFailures = 0; 1463d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 147184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 148184914b5SBarry Smith 1495e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1505e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15139e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1525e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1535e42470aSBarry Smith G = snes->work[1]; 1545e42470aSBarry Smith W = snes->work[2]; 1555e76c431SBarry Smith 1564c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1574c49b128SBarry Smith snes->iter = 0; 15875567043SBarry Smith snes->norm = 0.0; 1594c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16019717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1614936397dSBarry Smith if (snes->domainerror) { 1624936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1634936397dSBarry Smith PetscFunctionReturn(0); 1644936397dSBarry Smith } 1652613ca53SBarry Smith ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1662613ca53SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1672613ca53SBarry Smith ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 1682613ca53SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 16962cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1700f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1715e42470aSBarry Smith snes->norm = fnorm; 1720f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 17384c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 1747a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1753f149594SLisandro Dalcin 1763f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1773f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 17885385478SLisandro Dalcin /* test convergence */ 17985385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 18006ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 181d034289bSBarry Smith 1825e76c431SBarry Smith for (i=0; i<maxits; i++) { 1835e76c431SBarry Smith 184329e7e40SMatthew Knepley /* Call general purpose update function */ 185e7788613SBarry Smith if (snes->ops->update) { 186e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 187de8cb200SMatthew Knepley } 188329e7e40SMatthew Knepley 189ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1905334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 19194b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 19271f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 1933d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1943d4c4710SBarry Smith if (kspreason < 0) { 1953d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 196ef998cc9SBarry 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); 1973d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 198ab81109fSSatish Balay break; 1993d4c4710SBarry Smith } 2003d4c4710SBarry Smith } 2013d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2023f149594SLisandro Dalcin snes->linear_its += lits; 2033f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 20474637425SBarry Smith 2059c3ca13aSBarry Smith if (neP->precheckstep) { 206ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 2079c3ca13aSBarry Smith ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 2089c3ca13aSBarry Smith } 2099c3ca13aSBarry Smith 210b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 2111e2582c4SBarry Smith ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 21285471664SBarry Smith } 21374637425SBarry Smith 214ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 215ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 216e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 217ea4d3ed3SLois Curfman McInnes */ 21885385478SLisandro Dalcin ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 2193f149594SLisandro Dalcin ynorm = 1; gnorm = fnorm; 220dc357ad2SBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 2218f1a2a5eSBarry 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); 2224a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2234936397dSBarry Smith if (snes->domainerror) { 2244936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2254936397dSBarry Smith PetscFunctionReturn(0); 2264936397dSBarry Smith } 227a7cc72afSBarry Smith if (!lssucceed) { 22850ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 229ace3abfcSBarry Smith PetscBool ismin; 230647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2311e2582c4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 2328a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2333505fcc1SBarry Smith break; 2343505fcc1SBarry Smith } 23550ffb88aSMatthew Knepley } 23685385478SLisandro Dalcin /* Update function and solution vectors */ 23785385478SLisandro Dalcin fnorm = gnorm; 23885385478SLisandro Dalcin ierr = VecCopy(G,F);CHKERRQ(ierr); 23985385478SLisandro Dalcin ierr = VecCopy(W,X);CHKERRQ(ierr); 24085385478SLisandro Dalcin /* Monitor convergence */ 24185385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 24285385478SLisandro Dalcin snes->iter = i+1; 24385385478SLisandro Dalcin snes->norm = fnorm; 24485385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24585385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 2467a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 24785385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 24885385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 249e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2503f149594SLisandro Dalcin if (snes->reason) break; 25116c95cacSBarry Smith } 25252392280SLois Curfman McInnes if (i == maxits) { 253ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25485385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25552392280SLois Curfman McInnes } 2563a40ed3dSBarry Smith PetscFunctionReturn(0); 2575e76c431SBarry Smith } 25804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25904d965bbSLois Curfman McInnes /* 2604b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2614b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26204d965bbSLois Curfman McInnes 26304d965bbSLois Curfman McInnes Input Parameter: 26404d965bbSLois Curfman McInnes . snes - the SNES context 26504d965bbSLois Curfman McInnes . x - the solution vector 26604d965bbSLois Curfman McInnes 26704d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26804d965bbSLois Curfman McInnes 26904d965bbSLois Curfman McInnes Notes: 2704b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 27104d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27204d965bbSLois Curfman McInnes the call to SNESSolve(). 27304d965bbSLois Curfman McInnes */ 2744a2ae208SSatish Balay #undef __FUNCT__ 2754b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 276dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2775e76c431SBarry Smith { 278dfbe8321SBarry Smith PetscErrorCode ierr; 2793a40ed3dSBarry Smith 2803a40ed3dSBarry Smith PetscFunctionBegin; 28158c9b817SLisandro Dalcin ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 2823a40ed3dSBarry Smith PetscFunctionReturn(0); 2835e76c431SBarry Smith } 28404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2856b8b9a38SLisandro Dalcin 2866b8b9a38SLisandro Dalcin #undef __FUNCT__ 2876b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS" 2886b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes) 2896b8b9a38SLisandro Dalcin { 2906b8b9a38SLisandro Dalcin PetscErrorCode ierr; 2916b8b9a38SLisandro Dalcin 2926b8b9a38SLisandro Dalcin PetscFunctionBegin; 2936b8b9a38SLisandro Dalcin if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2946b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2956b8b9a38SLisandro Dalcin } 2966b8b9a38SLisandro Dalcin 29704d965bbSLois Curfman McInnes /* 2984b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2994b27c08aSLois Curfman McInnes with SNESCreate_LS(). 30004d965bbSLois Curfman McInnes 30104d965bbSLois Curfman McInnes Input Parameter: 30204d965bbSLois Curfman McInnes . snes - the SNES context 30304d965bbSLois Curfman McInnes 304de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30504d965bbSLois Curfman McInnes */ 3064a2ae208SSatish Balay #undef __FUNCT__ 3074b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 308dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3095e76c431SBarry Smith { 31094298653SBarry Smith SNES_LS *ls = (SNES_LS*) snes->data; 311dfbe8321SBarry Smith PetscErrorCode ierr; 3123a40ed3dSBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 3146b8b9a38SLisandro Dalcin ierr = SNESReset_LS(snes);CHKERRQ(ierr); 315649052a6SBarry Smith ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr); 3165c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 317557d3f75SLisandro Dalcin 318557d3f75SLisandro Dalcin /* clear composed functions */ 31958c9b817SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 32058c9b817SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);CHKERRQ(ierr); 321557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 322557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 323557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 324557d3f75SLisandro Dalcin 3253a40ed3dSBarry Smith PetscFunctionReturn(0); 3265e76c431SBarry Smith } 32704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3284a2ae208SSatish Balay #undef __FUNCT__ 32917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 33082bf6240SBarry Smith 3314b828684SBarry Smith /*@C 33217bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3335e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3345e76c431SBarry Smith to serve as a template and is not recommended for general use. 3355e76c431SBarry Smith 3363f9fe445SBarry Smith Logically Collective on SNES and Vec 337c7afd0dbSLois Curfman McInnes 3385e76c431SBarry Smith Input Parameters: 339c7afd0dbSLois Curfman McInnes + snes - nonlinear context 340af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3415e76c431SBarry Smith . x - current iterate 3425e76c431SBarry Smith . f - residual evaluated at x 3433c632250SBarry Smith . y - search direction 344dc357ad2SBarry Smith . fnorm - 2-norm of f 345dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 3465e76c431SBarry Smith 347c4a48953SLois Curfman McInnes Output Parameters: 348c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3493c632250SBarry Smith . w - new iterate 3505e76c431SBarry Smith . gnorm - 2-norm of g 3515e76c431SBarry Smith . ynorm - 2-norm of search length 352a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 353fee21e36SBarry Smith 354c4a48953SLois Curfman McInnes Options Database Key: 35517bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 356c4a48953SLois Curfman McInnes 35736851e7fSLois Curfman McInnes Level: advanced 35836851e7fSLois Curfman McInnes 35928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 36028ae5a21SLois Curfman McInnes 36117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 36217bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3635e76c431SBarry Smith @*/ 3647087cfbeSBarry 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) 3655e76c431SBarry Smith { 366dfbe8321SBarry Smith PetscErrorCode ierr; 3674b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 368ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3695334005bSBarry Smith 3703a40ed3dSBarry Smith PetscFunctionBegin; 371a7cc72afSBarry Smith *flag = PETSC_TRUE; 372d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 373a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 37479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3753c632250SBarry Smith if (neP->postcheckstep) { 3763c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3778f99978dSLois Curfman McInnes } 3783c632250SBarry Smith if (changed_y) { 37979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3803c632250SBarry Smith } 3814936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 3824936397dSBarry Smith if (!snes->domainerror) { 383a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 38462cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3854936397dSBarry Smith } 386d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3885e76c431SBarry Smith } 38904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3904a2ae208SSatish Balay #undef __FUNCT__ 39117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 39282bf6240SBarry Smith 39329e0b56fSBarry Smith /*@C 39417bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 39529e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 39629e0b56fSBarry Smith even compute the norm of the function or search direction; this 39729e0b56fSBarry Smith is intended only when you know the full step is fine and are 39829e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 399c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 400c7afd0dbSLois Curfman McInnes 4013f9fe445SBarry Smith Logically Collective on SNES and Vec 40229e0b56fSBarry Smith 40329e0b56fSBarry Smith Input Parameters: 404c7afd0dbSLois Curfman McInnes + snes - nonlinear context 405af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 40629e0b56fSBarry Smith . x - current iterate 40729e0b56fSBarry Smith . f - residual evaluated at x 4083c632250SBarry Smith . y - search direction 40929e0b56fSBarry Smith . w - work vector 410dc357ad2SBarry Smith . fnorm - 2-norm of f 411dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 41229e0b56fSBarry Smith 41329e0b56fSBarry Smith Output Parameters: 414c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4153c632250SBarry Smith . w - new iterate 41629e0b56fSBarry Smith . gnorm - not changed 41729e0b56fSBarry Smith . ynorm - not changed 418a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 419fee21e36SBarry Smith 42029e0b56fSBarry Smith Options Database Key: 42117bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 42229e0b56fSBarry Smith 4238cbba510SBarry Smith Notes: 42417bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 425ea56f5baSLois Curfman McInnes either the options 426ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 427ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 428329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 429329f5518SBarry Smith otherwise, the SNES solver will generate an error. 430329f5518SBarry Smith 431329f5518SBarry Smith During the final iteration this will not evaluate the function at 432329f5518SBarry Smith the solution point. This is to save a function evaluation while 433329f5518SBarry Smith using pseudo-timestepping. 4348cbba510SBarry Smith 435ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 436a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 437ea56f5baSLois Curfman McInnes correct, since they are not computed. 438ea56f5baSLois Curfman McInnes 439ea56f5baSLois Curfman McInnes Level: advanced 4408cbba510SBarry Smith 44129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 44229e0b56fSBarry Smith 44317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 44417bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 44529e0b56fSBarry Smith @*/ 4467087cfbeSBarry 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) 44729e0b56fSBarry Smith { 448dfbe8321SBarry Smith PetscErrorCode ierr; 4494b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 450ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 45129e0b56fSBarry Smith 4523a40ed3dSBarry Smith PetscFunctionBegin; 453a7cc72afSBarry Smith *flag = PETSC_TRUE; 454d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 45579f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4563c632250SBarry Smith if (neP->postcheckstep) { 4573c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4583c632250SBarry Smith } 4593c632250SBarry Smith if (changed_y) { 46079f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4618f99978dSLois Curfman McInnes } 462329f5518SBarry Smith 463329f5518SBarry Smith /* don't evaluate function the last time through */ 464329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4654936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 466329f5518SBarry Smith } 467d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 46929e0b56fSBarry Smith } 47004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4714a2ae208SSatish Balay #undef __FUNCT__ 47217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4734b828684SBarry Smith /*@C 47417bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4755e76c431SBarry Smith 476c7afd0dbSLois Curfman McInnes Collective on SNES 477c7afd0dbSLois Curfman McInnes 4785e76c431SBarry Smith Input Parameters: 479c7afd0dbSLois Curfman McInnes + snes - nonlinear context 480af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4815e76c431SBarry Smith . x - current iterate 4825e76c431SBarry Smith . f - residual evaluated at x 4833c632250SBarry Smith . y - search direction 4845e76c431SBarry Smith . w - work vector 485dc357ad2SBarry Smith . fnorm - 2-norm of f 486dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 4875e76c431SBarry Smith 488393d2d9aSLois Curfman McInnes Output Parameters: 489c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4903c632250SBarry Smith . w - new iterate 4915e76c431SBarry Smith . gnorm - 2-norm of g 4925e76c431SBarry Smith . ynorm - 2-norm of search length 493a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 494fee21e36SBarry Smith 495c4a48953SLois Curfman McInnes Options Database Key: 496e106eecfSBarry Smith + -snes_ls cubic - Activates SNESLineSearchCubic() 497e106eecfSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 498e106eecfSBarry 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) 499e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 500e106eecfSBarry Smith 501c4a48953SLois Curfman McInnes 5025e76c431SBarry Smith Notes: 5035e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 5045e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 5055e76c431SBarry Smith 50636851e7fSLois Curfman McInnes Level: advanced 50736851e7fSLois Curfman McInnes 50828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 50928ae5a21SLois Curfman McInnes 51017bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 5115e76c431SBarry Smith @*/ 5127087cfbeSBarry 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) 5135e76c431SBarry Smith { 514406556e6SLois Curfman McInnes /* 515406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 516406556e6SLois Curfman McInnes minimization problem: 517406556e6SLois Curfman McInnes min z(x): R^n -> R, 518406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 519406556e6SLois Curfman McInnes */ 520406556e6SLois Curfman McInnes 521e106eecfSBarry Smith PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 522e106eecfSBarry Smith PetscReal minlambda,lambda,lambdatemp; 523aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 524f5b6597dSBarry Smith PetscScalar cinitslope; 5256b5873e3SBarry Smith #endif 526dfbe8321SBarry Smith PetscErrorCode ierr; 527a7cc72afSBarry Smith PetscInt count; 5284b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 529ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 53094298653SBarry Smith MPI_Comm comm; 5315e76c431SBarry Smith 5323a40ed3dSBarry Smith PetscFunctionBegin; 53394298653SBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 534d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 535a7cc72afSBarry Smith *flag = PETSC_TRUE; 5365e76c431SBarry Smith 537cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 53875567043SBarry Smith if (*ynorm == 0.0) { 53994298653SBarry Smith if (neP->monitor) { 540649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 541649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 542649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 54394298653SBarry Smith } 544a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5453c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 546a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 547a7cc72afSBarry Smith *flag = PETSC_FALSE; 548ad922ac9SBarry Smith goto theend1; 54994a9d846SBarry Smith } 550e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 55194298653SBarry Smith if (neP->monitor) { 552649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 5538f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(neP->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr); 554649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 55594298653SBarry Smith } 556e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 557e106eecfSBarry Smith *ynorm = neP->maxstep; 5585e76c431SBarry Smith } 559ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 560e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 561a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 562aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 563a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 564329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5655e42470aSBarry Smith #else 566a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5675e42470aSBarry Smith #endif 5685e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5695e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5705e76c431SBarry Smith 571e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 57243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 573ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 57443e71028SBarry Smith *flag = PETSC_FALSE; 57543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 57643e71028SBarry Smith goto theend1; 57743e71028SBarry Smith } 5784936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 5794936397dSBarry Smith if (snes->domainerror) { 5804936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 5814936397dSBarry Smith PetscFunctionReturn(0); 58219717074SBarry Smith } 583313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 58462cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5858f1a2a5eSBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr); 586e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 58794298653SBarry Smith if (neP->monitor) { 588649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 5898f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr); 590649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 59194298653SBarry Smith } 592ad922ac9SBarry Smith goto theend1; 5935e76c431SBarry Smith } 5945e76c431SBarry Smith 5955e76c431SBarry Smith /* Fit points with quadratic */ 596313b5538SBarry Smith lambda = 1.0; 597a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5985e76c431SBarry Smith lambdaprev = lambda; 5995e76c431SBarry Smith gnormprev = *gnorm; 600329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 601ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 602ddd12767SBarry Smith else lambda = lambdatemp; 603275d25c3SBarry Smith 604e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 60543e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 606ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 60743e71028SBarry Smith *flag = PETSC_FALSE; 60843e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 60943e71028SBarry Smith goto theend1; 61043e71028SBarry Smith } 611f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6124936397dSBarry Smith if (snes->domainerror) { 6134936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6144936397dSBarry Smith PetscFunctionReturn(0); 61519717074SBarry Smith } 616cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 61762cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 61894298653SBarry Smith if (neP->monitor) { 619649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 6208f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr); 621649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 62294298653SBarry Smith } 623e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 62494298653SBarry Smith if (neP->monitor) { 625649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 626649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 627649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 62894298653SBarry Smith } 629ad922ac9SBarry Smith goto theend1; 6305e76c431SBarry Smith } 6315e76c431SBarry Smith 6325e76c431SBarry Smith /* Fit points with cubic */ 6335e76c431SBarry Smith count = 1; 634bb9195b6SBarry Smith while (PETSC_TRUE) { 635dc357ad2SBarry Smith if (lambda <= minlambda) { 63694298653SBarry Smith if (neP->monitor) { 637649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 638649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 639649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(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); 640649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 64194298653SBarry Smith } 64243e71028SBarry Smith *flag = PETSC_FALSE; 64343e71028SBarry Smith break; 6445e76c431SBarry Smith } 6455d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6465d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 647ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6482b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6495e76c431SBarry Smith d = b*b - 3*a*initslope; 6505e76c431SBarry Smith if (d < 0.0) d = 0.0; 6515e76c431SBarry Smith if (a == 0.0) { 6525e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6535e76c431SBarry Smith } else { 6545e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6555e76c431SBarry Smith } 6565e76c431SBarry Smith lambdaprev = lambda; 6575e76c431SBarry Smith gnormprev = *gnorm; 658329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 659329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6605e76c431SBarry Smith else lambda = lambdatemp; 661e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 66243e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 663ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 664ae15b995SBarry 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); 665ed98166eSMatthew Knepley *flag = PETSC_FALSE; 66643e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 667ed98166eSMatthew Knepley break; 668ed98166eSMatthew Knepley } 669f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6704936397dSBarry Smith if (snes->domainerror) { 6714936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6724936397dSBarry Smith PetscFunctionReturn(0); 67319717074SBarry Smith } 674cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 67562cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 676e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 67794298653SBarry Smith if (neP->monitor) { 6788f1a2a5eSBarry Smith ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr); 67994298653SBarry Smith } 68043e71028SBarry Smith break; 6812b022350SLois Curfman McInnes } else { 68294298653SBarry Smith if (neP->monitor) { 6838f1a2a5eSBarry 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); 68494298653SBarry Smith } 6855e76c431SBarry Smith } 6865e76c431SBarry Smith count++; 6875e76c431SBarry Smith } 688ad922ac9SBarry Smith theend1: 6898f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6903c632250SBarry Smith if (neP->postcheckstep && *flag) { 6913c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6923c632250SBarry Smith if (changed_y) { 69379f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6943c632250SBarry Smith } 6953c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 696f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 6974936397dSBarry Smith if (snes->domainerror) { 6984936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6994936397dSBarry Smith PetscFunctionReturn(0); 70019717074SBarry Smith } 7018f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 70262cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 703a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7048f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 705a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7068f99978dSLois Curfman McInnes } 7078f99978dSLois Curfman McInnes } 708d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 7093a40ed3dSBarry Smith PetscFunctionReturn(0); 7105e76c431SBarry Smith } 71104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7124a2ae208SSatish Balay #undef __FUNCT__ 71317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 7144b828684SBarry Smith /*@C 71517bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 7165e76c431SBarry Smith 717c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 718c7afd0dbSLois Curfman McInnes 7195e42470aSBarry Smith Input Parameters: 720c7afd0dbSLois Curfman McInnes + snes - the SNES context 721af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 7225e42470aSBarry Smith . x - current iterate 7235e42470aSBarry Smith . f - residual evaluated at x 7243c632250SBarry Smith . y - search direction 7255e42470aSBarry Smith . w - work vector 726dc357ad2SBarry Smith . fnorm - 2-norm of f 727dc357ad2SBarry Smith - xnorm - norm of x if known, otherwise 0 7285e42470aSBarry Smith 729c4a48953SLois Curfman McInnes Output Parameters: 7307f3332b4SBarry Smith + g - residual evaluated at new iterate w 731e106eecfSBarry Smith . w - new iterate (x + lambda*y) 7325e42470aSBarry Smith . gnorm - 2-norm of g 7335e42470aSBarry Smith . ynorm - 2-norm of search length 734a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 735fee21e36SBarry Smith 736ce9499c7SBarry Smith Options Database Keys: 737ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 738ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 739e106eecfSBarry 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) 740e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 741e106eecfSBarry Smith 7425e42470aSBarry Smith Notes: 74317bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7445e42470aSBarry Smith 74536851e7fSLois Curfman McInnes Level: advanced 74636851e7fSLois Curfman McInnes 747f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 748f59ffdedSLois Curfman McInnes 74917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7505e42470aSBarry Smith @*/ 7517087cfbeSBarry 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) 7525e76c431SBarry Smith { 753406556e6SLois Curfman McInnes /* 754406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 755406556e6SLois Curfman McInnes minimization problem: 756406556e6SLois Curfman McInnes min z(x): R^n -> R, 757406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 758406556e6SLois Curfman McInnes */ 759e106eecfSBarry Smith PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 760aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 761f5b6597dSBarry Smith PetscScalar cinitslope; 7626b5873e3SBarry Smith #endif 763dfbe8321SBarry Smith PetscErrorCode ierr; 764a7cc72afSBarry Smith PetscInt count; 7654b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 766ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7675e76c431SBarry Smith 7683a40ed3dSBarry Smith PetscFunctionBegin; 769d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 770a7cc72afSBarry Smith *flag = PETSC_TRUE; 7715e76c431SBarry Smith 7723505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 77363b9fa5eSBarry Smith if (*ynorm == 0.0) { 77494298653SBarry Smith if (neP->monitor) { 775649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 776649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 777649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 77894298653SBarry Smith } 779b37302c6SLois Curfman McInnes *gnorm = fnorm; 780e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 781b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 782a7cc72afSBarry Smith *flag = PETSC_FALSE; 783ad922ac9SBarry Smith goto theend2; 78494a9d846SBarry Smith } 785e106eecfSBarry Smith if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 786e106eecfSBarry Smith ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 787e106eecfSBarry Smith *ynorm = neP->maxstep; 7885e76c431SBarry Smith } 789ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 790e106eecfSBarry Smith minlambda = neP->minlambda/rellength; 791a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 792aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 793a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 794329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7955e42470aSBarry Smith #else 796a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7975e42470aSBarry Smith #endif 7985e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7995e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 8008f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr); 8015e42470aSBarry Smith 802e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 80343e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 804ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 80543e71028SBarry Smith *flag = PETSC_FALSE; 80643e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 80743e71028SBarry Smith goto theend2; 80843e71028SBarry Smith } 8094936397dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8104936397dSBarry Smith if (snes->domainerror) { 8114936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8124936397dSBarry Smith PetscFunctionReturn(0); 81319717074SBarry Smith } 814cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 81562cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 816e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 81794298653SBarry Smith if (neP->monitor) { 818649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 819649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 820649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 82194298653SBarry Smith } 822ad922ac9SBarry Smith goto theend2; 8235e42470aSBarry Smith } 8245e42470aSBarry Smith 8255e42470aSBarry Smith /* Fit points with quadratic */ 826313b5538SBarry Smith lambda = 1.0; 8275e42470aSBarry Smith count = 1; 8285ca10a37SBarry Smith while (PETSC_TRUE) { 8295e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 83094298653SBarry Smith if (neP->monitor) { 831649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 832649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 8338f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(neP->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); 834649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 83594298653SBarry Smith } 836e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 83743e71028SBarry Smith *flag = PETSC_FALSE; 83843e71028SBarry Smith break; 8395e42470aSBarry Smith } 840a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 841329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 842329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 843329f5518SBarry Smith else lambda = lambdatemp; 844275d25c3SBarry Smith 845e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 84643e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 847ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 848d9822059SBarry 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); 849ed98166eSMatthew Knepley *flag = PETSC_FALSE; 85043e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 851ed98166eSMatthew Knepley break; 852ed98166eSMatthew Knepley } 853f5b6597dSBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 8544936397dSBarry Smith if (snes->domainerror) { 8554936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8564936397dSBarry Smith PetscFunctionReturn(0); 85719717074SBarry Smith } 858cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 85962cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 860e106eecfSBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 86194298653SBarry Smith if (neP->monitor) { 862649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 8638f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(neP->monitor," Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr); 864649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 86594298653SBarry Smith } 86606259719SBarry Smith break; 8675e42470aSBarry Smith } 8685e42470aSBarry Smith count++; 8695e42470aSBarry Smith } 870ad922ac9SBarry Smith theend2: 8718f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8723c632250SBarry Smith if (neP->postcheckstep) { 8733c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8743c632250SBarry Smith if (changed_y) { 87579f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8763c632250SBarry Smith } 8773c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8783c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 8794936397dSBarry Smith if (snes->domainerror) { 8804936397dSBarry Smith ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8814936397dSBarry Smith PetscFunctionReturn(0); 88219717074SBarry Smith } 8838f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 884a9f58f88SMatthew Knepley ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8858f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 886a9f58f88SMatthew Knepley ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 88762cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8888f99978dSLois Curfman McInnes } 8898f99978dSLois Curfman McInnes } 890d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8913a40ed3dSBarry Smith PetscFunctionReturn(0); 8925e76c431SBarry Smith } 8932343ba6eSBarry Smith 89404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8954a2ae208SSatish Balay #undef __FUNCT__ 89617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 897c9e6a524SLois Curfman McInnes /*@C 89817bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8994b27c08aSLois Curfman McInnes by the method SNESLS. 9005e76c431SBarry Smith 9015e76c431SBarry Smith Input Parameters: 902c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 903af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 904c7afd0dbSLois Curfman McInnes - func - pointer to int function 9055e76c431SBarry Smith 9063f9fe445SBarry Smith Logically Collective on SNES 907fee21e36SBarry Smith 908c4a48953SLois Curfman McInnes Available Routines: 90917bae607SBarry Smith + SNESLineSearchCubic() - default line search 91017bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 91117bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 91217bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 9135e76c431SBarry Smith 914c4a48953SLois Curfman McInnes Options Database Keys: 9154b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 9164b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 917e106eecfSBarry 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) 918e106eecfSBarry Smith - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 919c4a48953SLois Curfman McInnes 9205e76c431SBarry Smith Calling sequence of func: 921bd208895SLois Curfman McInnes .vb 922ace3abfcSBarry 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) 923bd208895SLois Curfman McInnes .ve 9245e76c431SBarry Smith 9255e76c431SBarry Smith Input parameters for func: 926c7afd0dbSLois Curfman McInnes + snes - nonlinear context 927af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 9285e76c431SBarry Smith . x - current iterate 9295e76c431SBarry Smith . f - residual evaluated at x 9303c632250SBarry Smith . y - search direction 931c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9325e76c431SBarry Smith 9335e76c431SBarry Smith Output parameters for func: 934c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9353c632250SBarry Smith . w - new iterate 9365e76c431SBarry Smith . gnorm - 2-norm of g 9375e76c431SBarry Smith . ynorm - 2-norm of search length 938a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 939f59ffdedSLois Curfman McInnes 94036851e7fSLois Curfman McInnes Level: advanced 94136851e7fSLois Curfman McInnes 942f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 943f59ffdedSLois Curfman McInnes 94417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 94517bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9465e76c431SBarry Smith @*/ 9477087cfbeSBarry Smith PetscErrorCode SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx) 9485e76c431SBarry Smith { 9494ac538c5SBarry Smith PetscErrorCode ierr; 95082bf6240SBarry Smith 9513a40ed3dSBarry Smith PetscFunctionBegin; 9524ac538c5SBarry 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); 9533a40ed3dSBarry Smith PetscFunctionReturn(0); 9545e76c431SBarry Smith } 9558e019c35SBarry Smith 956ace3abfcSBarry 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*/ 95704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 958fb2e594dSBarry Smith EXTERN_C_BEGIN 9594a2ae208SSatish Balay #undef __FUNCT__ 96017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 9617087cfbeSBarry Smith PetscErrorCode SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 96282bf6240SBarry Smith { 96382bf6240SBarry Smith PetscFunctionBegin; 9644b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9654b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 96682bf6240SBarry Smith PetscFunctionReturn(0); 96782bf6240SBarry Smith } 968fb2e594dSBarry Smith EXTERN_C_END 96904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9704a2ae208SSatish Balay #undef __FUNCT__ 9713c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 972c8dd0c0dSLois Curfman McInnes /*@C 9733c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9744b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 975c8dd0c0dSLois Curfman McInnes 976c8dd0c0dSLois Curfman McInnes Input Parameters: 977c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9783c632250SBarry Smith . func - pointer to function 979c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 980c8dd0c0dSLois Curfman McInnes 9813f9fe445SBarry Smith Logically Collective on SNES 982c8dd0c0dSLois Curfman McInnes 983c8dd0c0dSLois Curfman McInnes Calling sequence of func: 984c8dd0c0dSLois Curfman McInnes .vb 985ace3abfcSBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool *changed_y,PetscBool *changed_w) 986c8dd0c0dSLois Curfman McInnes .ve 987b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 988b1ae27eaSLois Curfman McInnes on failure. 989c8dd0c0dSLois Curfman McInnes 990c8dd0c0dSLois Curfman McInnes Input parameters for func: 991c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 992c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9933c632250SBarry Smith . x - previous iterate 9943c632250SBarry Smith . y - new search direction and length 9953c632250SBarry Smith - w - current candidate iterate 996c8dd0c0dSLois Curfman McInnes 997c8dd0c0dSLois Curfman McInnes Output parameters for func: 9983c632250SBarry Smith + y - search direction (possibly changed) 9993c632250SBarry Smith . w - current iterate (possibly modified) 10003c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 10013c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 1002c8dd0c0dSLois Curfman McInnes 1003c8dd0c0dSLois Curfman McInnes Level: advanced 1004c8dd0c0dSLois Curfman McInnes 10059e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10069e247f21SBarry Smith 10073c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 10083c632250SBarry Smith 1009481b4421SBarry Smith On input w = x - y 10103c632250SBarry Smith 101117bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 1012b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 1013ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 10148f99978dSLois Curfman McInnes 101517bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 1016ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 1017ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 1018ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 10199e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 1020b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 10218f99978dSLois Curfman McInnes 1022c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 1023c8dd0c0dSLois Curfman McInnes 102417bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1025c8dd0c0dSLois Curfman McInnes @*/ 10267087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx) 1027c8dd0c0dSLois Curfman McInnes { 10284ac538c5SBarry Smith PetscErrorCode ierr; 1029c8dd0c0dSLois Curfman McInnes 1030c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10314ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 1032c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1033c8dd0c0dSLois Curfman McInnes } 103494298653SBarry Smith 10359c3ca13aSBarry Smith #undef __FUNCT__ 10369c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10379c3ca13aSBarry Smith /*@C 10389c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10397e4bb74cSBarry Smith before the line search is called. 10409c3ca13aSBarry Smith 10419c3ca13aSBarry Smith Input Parameters: 10429c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10439c3ca13aSBarry Smith . func - pointer to function 10449c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10459c3ca13aSBarry Smith 10463f9fe445SBarry Smith Logically Collective on SNES 10479c3ca13aSBarry Smith 10489c3ca13aSBarry Smith Calling sequence of func: 10499c3ca13aSBarry Smith .vb 1050ace3abfcSBarry Smith int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool *changed_y) 10519c3ca13aSBarry Smith .ve 10529c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10539c3ca13aSBarry Smith on failure. 10549c3ca13aSBarry Smith 10559c3ca13aSBarry Smith Input parameters for func: 10569c3ca13aSBarry Smith + snes - nonlinear context 10579c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10589c3ca13aSBarry Smith . x - previous iterate 10599c3ca13aSBarry Smith - y - new search direction and length 10609c3ca13aSBarry Smith 10619c3ca13aSBarry Smith Output parameters for func: 10629c3ca13aSBarry Smith + y - search direction (possibly changed) 10639c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10649c3ca13aSBarry Smith 10659c3ca13aSBarry Smith Level: advanced 10669c3ca13aSBarry Smith 10679c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10689c3ca13aSBarry Smith 10699c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10709c3ca13aSBarry Smith 10717e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 10729c3ca13aSBarry Smith @*/ 10737087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx) 10749c3ca13aSBarry Smith { 10754ac538c5SBarry Smith PetscErrorCode ierr; 10769c3ca13aSBarry Smith 10779c3ca13aSBarry Smith PetscFunctionBegin; 10784ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 10799c3ca13aSBarry Smith PetscFunctionReturn(0); 10809c3ca13aSBarry Smith } 10819c3ca13aSBarry Smith 108294298653SBarry Smith #undef __FUNCT__ 108394298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor" 108494298653SBarry Smith /*@C 108594298653SBarry Smith SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search 108694298653SBarry Smith 108794298653SBarry Smith Input Parameters: 108894298653SBarry Smith + snes - nonlinear context obtained from SNESCreate() 108994298653SBarry Smith - flg - PETSC_TRUE to monitor the line search 109094298653SBarry Smith 10913f9fe445SBarry Smith Logically Collective on SNES 109294298653SBarry Smith 109394298653SBarry Smith Options Database: 109494298653SBarry Smith . -snes_ls_monitor 109594298653SBarry Smith 109694298653SBarry Smith Level: intermediate 109794298653SBarry Smith 109894298653SBarry Smith 109994298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 110094298653SBarry Smith @*/ 11017087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetMonitor(SNES snes,PetscBool flg) 110294298653SBarry Smith { 11034ac538c5SBarry Smith PetscErrorCode ierr; 110494298653SBarry Smith 110594298653SBarry Smith PetscFunctionBegin; 11064ac538c5SBarry Smith ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr); 110794298653SBarry Smith PetscFunctionReturn(0); 110894298653SBarry Smith } 110994298653SBarry Smith 1110c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1111ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/ 1112ace3abfcSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *); /* force argument to next function to not be extern C*/ 1113c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 11144a2ae208SSatish Balay #undef __FUNCT__ 11153c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 11167087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1117c8dd0c0dSLois Curfman McInnes { 1118c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 11193c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 11203c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1121c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1122c8dd0c0dSLois Curfman McInnes } 1123c8dd0c0dSLois Curfman McInnes EXTERN_C_END 11249c3ca13aSBarry Smith 11259c3ca13aSBarry Smith EXTERN_C_BEGIN 11269c3ca13aSBarry Smith #undef __FUNCT__ 11279c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 11287087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 11299c3ca13aSBarry Smith { 11309c3ca13aSBarry Smith PetscFunctionBegin; 11319c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 11329c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 11339c3ca13aSBarry Smith PetscFunctionReturn(0); 11349c3ca13aSBarry Smith } 11359c3ca13aSBarry Smith EXTERN_C_END 1136329e7e40SMatthew Knepley 113794298653SBarry Smith EXTERN_C_BEGIN 113894298653SBarry Smith #undef __FUNCT__ 113994298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS" 11407087cfbeSBarry Smith PetscErrorCode SNESLineSearchSetMonitor_LS(SNES snes,PetscBool flg) 114194298653SBarry Smith { 114294298653SBarry Smith SNES_LS *ls = (SNES_LS*)snes->data; 114394298653SBarry Smith PetscErrorCode ierr; 114494298653SBarry Smith 114594298653SBarry Smith PetscFunctionBegin; 114694298653SBarry Smith if (flg && !ls->monitor) { 1147649052a6SBarry Smith ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&ls->monitor);CHKERRQ(ierr); 114894298653SBarry Smith } else if (!flg && ls->monitor) { 1149649052a6SBarry Smith ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr); 115094298653SBarry Smith } 115194298653SBarry Smith PetscFunctionReturn(0); 115294298653SBarry Smith } 115394298653SBarry Smith EXTERN_C_END 115494298653SBarry Smith 1155329e7e40SMatthew Knepley /* 11564b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 115704d965bbSLois Curfman McInnes 115804d965bbSLois Curfman McInnes Input Parameters: 115904d965bbSLois Curfman McInnes . SNES - the SNES context 116004d965bbSLois Curfman McInnes . viewer - visualization context 116104d965bbSLois Curfman McInnes 116204d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 116304d965bbSLois Curfman McInnes */ 11644a2ae208SSatish Balay #undef __FUNCT__ 11654b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11666849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1167a935fc98SLois Curfman McInnes { 11684b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11692fc52814SBarry Smith const char *cstr; 1170dfbe8321SBarry Smith PetscErrorCode ierr; 1171ace3abfcSBarry Smith PetscBool iascii; 1172a935fc98SLois Curfman McInnes 11733a40ed3dSBarry Smith PetscFunctionBegin; 11742692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 117532077d6dSBarry Smith if (iascii) { 117617bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 117717bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 117817bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 117919bcc07fSBarry Smith else cstr = "unknown"; 1180b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 11818f1a2a5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)ls->alpha,(double)ls->maxstep,(double)ls->minlambda);CHKERRQ(ierr); 118219bcc07fSBarry Smith } 11833a40ed3dSBarry Smith PetscFunctionReturn(0); 1184a935fc98SLois Curfman McInnes } 118504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 118604d965bbSLois Curfman McInnes /* 11874b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 118804d965bbSLois Curfman McInnes 118904d965bbSLois Curfman McInnes Input Parameter: 119004d965bbSLois Curfman McInnes . snes - the SNES context 119104d965bbSLois Curfman McInnes 119204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 119304d965bbSLois Curfman McInnes */ 11944a2ae208SSatish Balay #undef __FUNCT__ 11954b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11966849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11975e42470aSBarry Smith { 11984b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1199dfbe8321SBarry Smith PetscErrorCode ierr; 1200*490ca623SBarry Smith SNESLineSearchType indx; 1201ace3abfcSBarry Smith PetscBool flg,set; 12025e42470aSBarry Smith 12033a40ed3dSBarry Smith PetscFunctionBegin; 1204b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 12054b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 12064b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1207e106eecfSBarry Smith ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1208acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 120994298653SBarry Smith if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1210186905e3SBarry Smith 1211*490ca623SBarry Smith ierr = PetscOptionsEnum("-snes_ls","Line search used","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)SNES_LS_CUBIC,(PetscEnum*)&indx,&flg);CHKERRQ(ierr); 121225cdf11fSBarry Smith if (flg) { 121322e36657SBarry Smith switch (indx) { 1214*490ca623SBarry Smith case SNES_LS_BASIC: 121517bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1216b49fd9e1SBarry Smith break; 1217*490ca623SBarry Smith case SNES_LS_BASIC_NONORMS: 121817bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1219b49fd9e1SBarry Smith break; 1220*490ca623SBarry Smith case SNES_LS_QUADRATIC: 122117bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1222b49fd9e1SBarry Smith break; 1223*490ca623SBarry Smith case SNES_LS_CUBIC: 122417bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1225b49fd9e1SBarry Smith break; 12265e42470aSBarry Smith } 12275e42470aSBarry Smith } 1228b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 12293a40ed3dSBarry Smith PetscFunctionReturn(0); 12305e42470aSBarry Smith } 12314827ddcaSBarry Smith 12324827ddcaSBarry Smith EXTERN_C_BEGIN 12334827ddcaSBarry Smith extern PetscErrorCode SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal); 12344827ddcaSBarry Smith EXTERN_C_END 12354827ddcaSBarry Smith 123604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1237ebe3b25bSBarry Smith /*MC 1238ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 123904d965bbSLois Curfman McInnes 1240ebe3b25bSBarry Smith Options Database: 1241ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1242ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1243e106eecfSBarry 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) 1244acbee50cSBarry Smith . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1245acbee50cSBarry Smith - -snes_ls_monitor - print information about progress of line searches 1246acbee50cSBarry Smith 124704d965bbSLois Curfman McInnes 1248ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 124904d965bbSLois Curfman McInnes 1250ee3001cbSBarry Smith Level: beginner 1251ee3001cbSBarry Smith 125217bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 125317bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1254b3dd4ab5SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1255ebe3b25bSBarry Smith 1256ebe3b25bSBarry Smith M*/ 1257fb2e594dSBarry Smith EXTERN_C_BEGIN 12584a2ae208SSatish Balay #undef __FUNCT__ 12594b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 12607087cfbeSBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 12615e42470aSBarry Smith { 1262dfbe8321SBarry Smith PetscErrorCode ierr; 12634b27c08aSLois Curfman McInnes SNES_LS *neP; 12645e42470aSBarry Smith 12653a40ed3dSBarry Smith PetscFunctionBegin; 1266e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 1267e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 1268e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 1269e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 1270e7788613SBarry Smith snes->ops->view = SNESView_LS; 12716b8b9a38SLisandro Dalcin snes->ops->reset = SNESReset_LS; 12725e42470aSBarry Smith 127338f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 12745e42470aSBarry Smith snes->data = (void*)neP; 12755e42470aSBarry Smith neP->alpha = 1.e-4; 12765e42470aSBarry Smith neP->maxstep = 1.e8; 1277e106eecfSBarry Smith neP->minlambda = 1.e-12; 127817bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1279c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12803c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12813c632250SBarry Smith neP->postcheck = PETSC_NULL; 12823c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12833c632250SBarry Smith neP->precheck = PETSC_NULL; 128482bf6240SBarry Smith 128594298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr); 12864827ddcaSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr); 128794298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 128894298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 128994298653SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 129082bf6240SBarry Smith 12913a40ed3dSBarry Smith PetscFunctionReturn(0); 12925e42470aSBarry Smith } 1293fb2e594dSBarry Smith EXTERN_C_END 129482bf6240SBarry Smith 129582bf6240SBarry Smith 129682bf6240SBarry Smith 1297