163dd3a1aSKris Buschelman #define PETSCSNES_DLL 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 58a5d9ee4SBarry Smith /* 68a5d9ee4SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the function, 78a5d9ee4SBarry Smith but not a zero. In the case when one cannot compute J^T F we use the fact that 836669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 936669109SBarry Smith for this trick. 108a5d9ee4SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 16dfbe8321SBarry Smith PetscErrorCode ierr; 1736669109SBarry Smith PetscTruth hastranspose; 188a5d9ee4SBarry Smith 198a5d9ee4SBarry Smith PetscFunctionBegin; 208a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2136669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2236669109SBarry Smith if (hastranspose) { 238a5d9ee4SBarry Smith /* Compute || J^T F|| */ 248a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 258a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 26ae15b995SBarry Smith ierr = PetscInfo1(0,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 30ea709b57SSatish Balay PetscScalar result; 3136669109SBarry Smith PetscReal wnorm; 3236669109SBarry Smith 33abb0e124SMatthew Knepley ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDestroy(work);CHKERRQ(ierr); 3936669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 40ae15b995SBarry Smith ierr = PetscInfo1(0,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 475ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4874637425SBarry Smith */ 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5274637425SBarry Smith { 5374637425SBarry Smith PetscReal a1,a2; 54dfbe8321SBarry Smith PetscErrorCode ierr; 5574637425SBarry Smith PetscTruth hastranspose; 5674637425SBarry Smith 5774637425SBarry Smith PetscFunctionBegin; 5874637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 5974637425SBarry Smith if (hastranspose) { 6074637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6179f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6274637425SBarry Smith 6374637425SBarry Smith /* Compute || J^T W|| */ 6474637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6574637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 679e247f21SBarry Smith if (a1) { 68ae15b995SBarry Smith ierr = PetscInfo1(0,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 6985471664SBarry Smith } 7074637425SBarry Smith } 7174637425SBarry Smith PetscFunctionReturn(0); 7274637425SBarry Smith } 7374637425SBarry Smith 7404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7504d965bbSLois Curfman McInnes 7604d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7794b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7804d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7904d965bbSLois Curfman McInnes respectively. 8004d965bbSLois Curfman McInnes 81fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8204d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8304d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 84fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8504d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8604d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 874b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 884b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 8904d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9004d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9104d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9204d965bbSLois Curfman McInnes for all nonlinear solvers. 9304d965bbSLois Curfman McInnes 9404d965bbSLois Curfman McInnes Another key routine is: 9504d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9604d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9704d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9804d965bbSLois Curfman McInnes SNESSolve() if necessary. 9904d965bbSLois Curfman McInnes 10004d965bbSLois Curfman McInnes Additional basic routines are: 10104d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10204d965bbSLois Curfman McInnes have actually been used. 10304d965bbSLois Curfman McInnes These are called by application codes via the interface routines 104186905e3SBarry Smith SNESView(). 10504d965bbSLois Curfman McInnes 10604d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10704d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10804d965bbSLois Curfman McInnes above description applies to these categories also. 10904d965bbSLois Curfman McInnes 11004d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1115e42470aSBarry Smith /* 1124b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11304d965bbSLois Curfman McInnes method with a line search. 1145e76c431SBarry Smith 11504d965bbSLois Curfman McInnes Input Parameters: 11604d965bbSLois Curfman McInnes . snes - the SNES context 1175e76c431SBarry Smith 11804d965bbSLois Curfman McInnes Output Parameter: 119c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12004d965bbSLois Curfman McInnes 12104d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1225e76c431SBarry Smith 1235e76c431SBarry Smith Notes: 1245e76c431SBarry Smith This implements essentially a truncated Newton method with a 1255e76c431SBarry Smith line search. By default a cubic backtracking line search 1265e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1275e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 128393d2d9aSLois Curfman McInnes and Schnabel. 1295e76c431SBarry Smith */ 1304a2ae208SSatish Balay #undef __FUNCT__ 1314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 132dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1335e76c431SBarry Smith { 1344b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1356849ba73SBarry Smith PetscErrorCode ierr; 136a7cc72afSBarry Smith PetscInt maxits,i,lits; 137a7cc72afSBarry Smith PetscTruth lssucceed; 138112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1405e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 1413d4c4710SBarry Smith KSPConvergedReason kspreason; 1425e76c431SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 1443d4c4710SBarry Smith snes->numFailures = 0; 1453d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 146184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 147184914b5SBarry Smith 1485e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1495e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15039e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1515e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1525e42470aSBarry Smith G = snes->work[1]; 1535e42470aSBarry Smith W = snes->work[2]; 1545e76c431SBarry Smith 1554c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1564c49b128SBarry Smith snes->iter = 0; 1574c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 15819717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 159cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 16043e71028SBarry Smith if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1610f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1625e42470aSBarry Smith snes->norm = fnorm; 1630f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16484c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16594a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1663f149594SLisandro Dalcin 1673f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1683f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 1693f149594SLisandro Dalcin if (snes->ops->converged) { 170e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1713f149594SLisandro Dalcin } 17206ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 173d034289bSBarry Smith 1745e76c431SBarry Smith for (i=0; i<maxits; i++) { 1755e76c431SBarry Smith 176329e7e40SMatthew Knepley /* Call general purpose update function */ 177e7788613SBarry Smith if (snes->ops->update) { 178e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 179de8cb200SMatthew Knepley } 180329e7e40SMatthew Knepley 181ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1825334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 18394b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 18471f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 1853d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1863d4c4710SBarry Smith if (kspreason < 0) { 1873d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 188*ef998cc9SBarry 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); 1893d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 190ab81109fSSatish Balay break; 1913d4c4710SBarry Smith } 1923d4c4710SBarry Smith } 1933d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1943f149594SLisandro Dalcin snes->linear_its += lits; 1953f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 19674637425SBarry Smith 1979c3ca13aSBarry Smith if (neP->precheckstep) { 1989c3ca13aSBarry Smith PetscTruth changed_y = PETSC_FALSE; 1999c3ca13aSBarry Smith ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 2009c3ca13aSBarry Smith } 2019c3ca13aSBarry Smith 202b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 20374637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 20485471664SBarry Smith } 20574637425SBarry Smith 206ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 207ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 208e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 209ea4d3ed3SLois Curfman McInnes */ 21081b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 2113f149594SLisandro Dalcin ynorm = 1; gnorm = fnorm; 212a7cc72afSBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 213ae15b995SBarry Smith ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 214a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 2153c632250SBarry Smith TMP = X; X = W; snes->vec_sol_always = X; W = TMP; 216a3b891d8SBarry Smith fnorm = gnorm; 2174a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 218a3b891d8SBarry Smith 2195ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2205ed2d596SBarry Smith snes->iter = i+1; 2215ed2d596SBarry Smith snes->norm = fnorm; 2225ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2235ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2245ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2255ed2d596SBarry Smith 226a7cc72afSBarry Smith if (!lssucceed) { 2278a5d9ee4SBarry Smith PetscTruth ismin; 22850ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2293505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2308a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2318a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2323505fcc1SBarry Smith break; 2333505fcc1SBarry Smith } 23450ffb88aSMatthew Knepley } 2355e76c431SBarry Smith 2365e76c431SBarry Smith /* Test for convergence */ 237e7788613SBarry Smith if (snes->ops->converged) { 23829e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 239e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2403f149594SLisandro Dalcin if (snes->reason) break; 24116c95cacSBarry Smith } 24229e0b56fSBarry Smith } 24339e2f89bSBarry Smith if (X != snes->vec_sol) { 244393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 245e15f7bb5SBarry Smith } 246e15f7bb5SBarry Smith if (F != snes->vec_func) { 247e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 248e15f7bb5SBarry Smith } 24939e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 25039e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 25152392280SLois Curfman McInnes if (i == maxits) { 252ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 2533505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 25452392280SLois Curfman McInnes } 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 2565e76c431SBarry Smith } 25704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 25804d965bbSLois Curfman McInnes /* 2594b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2604b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26104d965bbSLois Curfman McInnes 26204d965bbSLois Curfman McInnes Input Parameter: 26304d965bbSLois Curfman McInnes . snes - the SNES context 26404d965bbSLois Curfman McInnes . x - the solution vector 26504d965bbSLois Curfman McInnes 26604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26704d965bbSLois Curfman McInnes 26804d965bbSLois Curfman McInnes Notes: 2694b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 27004d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27104d965bbSLois Curfman McInnes the call to SNESSolve(). 27204d965bbSLois Curfman McInnes */ 2734a2ae208SSatish Balay #undef __FUNCT__ 2744b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 275dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2765e76c431SBarry Smith { 277dfbe8321SBarry Smith PetscErrorCode ierr; 2783a40ed3dSBarry Smith 2793a40ed3dSBarry Smith PetscFunctionBegin; 280e74804ceSBarry Smith if (!snes->work) { 28181b6cf68SLois Curfman McInnes snes->nwork = 4; 28201bab653SMatthew Knepley ierr = VecDuplicateVecs(snes->vec_sol_always,snes->nwork,&snes->work);CHKERRQ(ierr); 283efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 28481b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 285e74804ceSBarry Smith } 2863a40ed3dSBarry Smith PetscFunctionReturn(0); 2875e76c431SBarry Smith } 28804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 28904d965bbSLois Curfman McInnes /* 2904b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2914b27c08aSLois Curfman McInnes with SNESCreate_LS(). 29204d965bbSLois Curfman McInnes 29304d965bbSLois Curfman McInnes Input Parameter: 29404d965bbSLois Curfman McInnes . snes - the SNES context 29504d965bbSLois Curfman McInnes 296de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 29704d965bbSLois Curfman McInnes */ 2984a2ae208SSatish Balay #undef __FUNCT__ 2994b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 300dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3015e76c431SBarry Smith { 302dfbe8321SBarry Smith PetscErrorCode ierr; 3033a40ed3dSBarry Smith 3043a40ed3dSBarry Smith PetscFunctionBegin; 3055baf8537SBarry Smith if (snes->nwork) { 3064b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 30721c89e3eSBarry Smith } 3085c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 309557d3f75SLisandro Dalcin 310557d3f75SLisandro Dalcin /* clear composed functions */ 311557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 312557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 313557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 314557d3f75SLisandro Dalcin 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 3165e76c431SBarry Smith } 31704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3184a2ae208SSatish Balay #undef __FUNCT__ 31917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo" 32082bf6240SBarry Smith 3214b828684SBarry Smith /*@C 32217bae607SBarry Smith SNESLineSearchNo - This routine is not a line search at all; 3235e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3245e76c431SBarry Smith to serve as a template and is not recommended for general use. 3255e76c431SBarry Smith 326c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 327c7afd0dbSLois Curfman McInnes 3285e76c431SBarry Smith Input Parameters: 329c7afd0dbSLois Curfman McInnes + snes - nonlinear context 330af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3315e76c431SBarry Smith . x - current iterate 3325e76c431SBarry Smith . f - residual evaluated at x 3333c632250SBarry Smith . y - search direction 3345e76c431SBarry Smith . w - work vector 335c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3365e76c431SBarry Smith 337c4a48953SLois Curfman McInnes Output Parameters: 338c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3393c632250SBarry Smith . w - new iterate 3405e76c431SBarry Smith . gnorm - 2-norm of g 3415e76c431SBarry Smith . ynorm - 2-norm of search length 342a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 343fee21e36SBarry Smith 344c4a48953SLois Curfman McInnes Options Database Key: 34517bae607SBarry Smith . -snes_ls basic - Activates SNESLineSearchNo() 346c4a48953SLois Curfman McInnes 34736851e7fSLois Curfman McInnes Level: advanced 34836851e7fSLois Curfman McInnes 34928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 35028ae5a21SLois Curfman McInnes 35117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 35217bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms() 3535e76c431SBarry Smith @*/ 35417bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 3555e76c431SBarry Smith { 356dfbe8321SBarry Smith PetscErrorCode ierr; 3574b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3583c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 3595334005bSBarry Smith 3603a40ed3dSBarry Smith PetscFunctionBegin; 361a7cc72afSBarry Smith *flag = PETSC_TRUE; 362d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 363a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 36479f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3653c632250SBarry Smith if (neP->postcheckstep) { 3663c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 3678f99978dSLois Curfman McInnes } 3683c632250SBarry Smith if (changed_y) { 36979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 3703c632250SBarry Smith } 3713c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 372d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 37319717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 37419717074SBarry Smith } 375d5e45103SBarry Smith CHKERRQ(ierr); 376d5e45103SBarry Smith 377a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 37843e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 379d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3803a40ed3dSBarry Smith PetscFunctionReturn(0); 3815e76c431SBarry Smith } 38204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3834a2ae208SSatish Balay #undef __FUNCT__ 38417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms" 38582bf6240SBarry Smith 38629e0b56fSBarry Smith /*@C 38717bae607SBarry Smith SNESLineSearchNoNorms - This routine is not a line search at 38829e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 38929e0b56fSBarry Smith even compute the norm of the function or search direction; this 39029e0b56fSBarry Smith is intended only when you know the full step is fine and are 39129e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 392c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 393c7afd0dbSLois Curfman McInnes 394c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 39529e0b56fSBarry Smith 39629e0b56fSBarry Smith Input Parameters: 397c7afd0dbSLois Curfman McInnes + snes - nonlinear context 398af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 39929e0b56fSBarry Smith . x - current iterate 40029e0b56fSBarry Smith . f - residual evaluated at x 4013c632250SBarry Smith . y - search direction 40229e0b56fSBarry Smith . w - work vector 403c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 40429e0b56fSBarry Smith 40529e0b56fSBarry Smith Output Parameters: 406c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4073c632250SBarry Smith . w - new iterate 40829e0b56fSBarry Smith . gnorm - not changed 40929e0b56fSBarry Smith . ynorm - not changed 410a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 411fee21e36SBarry Smith 41229e0b56fSBarry Smith Options Database Key: 41317bae607SBarry Smith . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 41429e0b56fSBarry Smith 4158cbba510SBarry Smith Notes: 41617bae607SBarry Smith SNESLineSearchNoNorms() must be used in conjunction with 417ea56f5baSLois Curfman McInnes either the options 418ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 419ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 420329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 421329f5518SBarry Smith otherwise, the SNES solver will generate an error. 422329f5518SBarry Smith 423329f5518SBarry Smith During the final iteration this will not evaluate the function at 424329f5518SBarry Smith the solution point. This is to save a function evaluation while 425329f5518SBarry Smith using pseudo-timestepping. 4268cbba510SBarry Smith 427ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 428a6570f20SBarry Smith SNESMonitorDefault() (as activated via -snes_monitor) will not be 429ea56f5baSLois Curfman McInnes correct, since they are not computed. 430ea56f5baSLois Curfman McInnes 431ea56f5baSLois Curfman McInnes Level: advanced 4328cbba510SBarry Smith 43329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 43429e0b56fSBarry Smith 43517bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 43617bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNo() 43729e0b56fSBarry Smith @*/ 43817bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 43929e0b56fSBarry Smith { 440dfbe8321SBarry Smith PetscErrorCode ierr; 4414b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4423c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 44329e0b56fSBarry Smith 4443a40ed3dSBarry Smith PetscFunctionBegin; 445a7cc72afSBarry Smith *flag = PETSC_TRUE; 446d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 44779f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4483c632250SBarry Smith if (neP->postcheckstep) { 4493c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 4503c632250SBarry Smith } 4513c632250SBarry Smith if (changed_y) { 45279f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 4538f99978dSLois Curfman McInnes } 454329f5518SBarry Smith 455329f5518SBarry Smith /* don't evaluate function the last time through */ 456329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 4573c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 458d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 45919717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 46019717074SBarry Smith } 461d5e45103SBarry Smith CHKERRQ(ierr); 462329f5518SBarry Smith } 463d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4643a40ed3dSBarry Smith PetscFunctionReturn(0); 46529e0b56fSBarry Smith } 46604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4674a2ae208SSatish Balay #undef __FUNCT__ 46817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic" 4694b828684SBarry Smith /*@C 47017bae607SBarry Smith SNESLineSearchCubic - Performs a cubic line search (default line search method). 4715e76c431SBarry Smith 472c7afd0dbSLois Curfman McInnes Collective on SNES 473c7afd0dbSLois Curfman McInnes 4745e76c431SBarry Smith Input Parameters: 475c7afd0dbSLois Curfman McInnes + snes - nonlinear context 476af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4775e76c431SBarry Smith . x - current iterate 4785e76c431SBarry Smith . f - residual evaluated at x 4793c632250SBarry Smith . y - search direction 4805e76c431SBarry Smith . w - work vector 481c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4825e76c431SBarry Smith 483393d2d9aSLois Curfman McInnes Output Parameters: 484c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4853c632250SBarry Smith . w - new iterate 4865e76c431SBarry Smith . gnorm - 2-norm of g 4875e76c431SBarry Smith . ynorm - 2-norm of search length 488a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 489fee21e36SBarry Smith 490c4a48953SLois Curfman McInnes Options Database Key: 49117bae607SBarry Smith $ -snes_ls cubic - Activates SNESLineSearchCubic() 492c4a48953SLois Curfman McInnes 4935e76c431SBarry Smith Notes: 4945e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4955e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4965e76c431SBarry Smith 49736851e7fSLois Curfman McInnes Level: advanced 49836851e7fSLois Curfman McInnes 49928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 50028ae5a21SLois Curfman McInnes 50117bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 5025e76c431SBarry Smith @*/ 50317bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 5045e76c431SBarry Smith { 505406556e6SLois Curfman McInnes /* 506406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 507406556e6SLois Curfman McInnes minimization problem: 508406556e6SLois Curfman McInnes min z(x): R^n -> R, 509406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 510406556e6SLois Curfman McInnes */ 511406556e6SLois Curfman McInnes 5125ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 513275d25c3SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp; 514aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 515ea709b57SSatish Balay PetscScalar cinitslope,clambda; 5166b5873e3SBarry Smith #endif 517dfbe8321SBarry Smith PetscErrorCode ierr; 518a7cc72afSBarry Smith PetscInt count; 5194b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 52079f36460SBarry Smith PetscScalar scale; 5213c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 5225e76c431SBarry Smith 5233a40ed3dSBarry Smith PetscFunctionBegin; 524d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 525a7cc72afSBarry Smith *flag = PETSC_TRUE; 5265e76c431SBarry Smith alpha = neP->alpha; 5275e76c431SBarry Smith maxstep = neP->maxstep; 5285e76c431SBarry Smith steptol = neP->steptol; 5295e76c431SBarry Smith 530cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 5319e247f21SBarry Smith if (!*ynorm) { 532ae15b995SBarry Smith ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 533a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 5343c632250SBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 535a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 536a7cc72afSBarry Smith *flag = PETSC_FALSE; 537ad922ac9SBarry Smith goto theend1; 53894a9d846SBarry Smith } 5395e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5405e42470aSBarry Smith scale = maxstep/(*ynorm); 541ae15b995SBarry Smith ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr); 5422dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 5435e76c431SBarry Smith *ynorm = maxstep; 5445e76c431SBarry Smith } 545ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5465ca10a37SBarry Smith minlambda = steptol/rellength; 547a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 548aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 549a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 550329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5515e42470aSBarry Smith #else 552a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5535e42470aSBarry Smith #endif 5545e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5555e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5565e76c431SBarry Smith 557e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 55843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 559ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 56043e71028SBarry Smith *flag = PETSC_FALSE; 56143e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 56243e71028SBarry Smith goto theend1; 56343e71028SBarry Smith } 56419717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 565d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 56619717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 56719717074SBarry Smith } 568d5e45103SBarry Smith CHKERRQ(ierr); 569313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 57043e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5715d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 572ae15b995SBarry Smith ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 573ad922ac9SBarry Smith goto theend1; 5745e76c431SBarry Smith } 5755e76c431SBarry Smith 5765e76c431SBarry Smith /* Fit points with quadratic */ 577313b5538SBarry Smith lambda = 1.0; 578a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5795e76c431SBarry Smith lambdaprev = lambda; 5805e76c431SBarry Smith gnormprev = *gnorm; 581329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 582ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 583ddd12767SBarry Smith else lambda = lambdatemp; 584275d25c3SBarry Smith 585aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 586a23f76dfSSatish Balay clambda = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 5875e42470aSBarry Smith #else 588e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 5895e42470aSBarry Smith #endif 59043e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 591ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 59243e71028SBarry Smith *flag = PETSC_FALSE; 59343e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 59443e71028SBarry Smith goto theend1; 59543e71028SBarry Smith } 59619717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 597d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 59819717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 59919717074SBarry Smith } 600d5e45103SBarry Smith CHKERRQ(ierr); 601cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 60243e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6035ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 604ae15b995SBarry Smith ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 605ad922ac9SBarry Smith goto theend1; 6065e76c431SBarry Smith } 6075e76c431SBarry Smith 6085e76c431SBarry Smith /* Fit points with cubic */ 6095e76c431SBarry Smith count = 1; 6108229bfc2SMatthew Knepley while (count < 10000) { 6115e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 612ae15b995SBarry Smith ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 613ae15b995SBarry 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); 61443e71028SBarry Smith *flag = PETSC_FALSE; 61543e71028SBarry Smith break; 6165e76c431SBarry Smith } 6175d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6185d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 619ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6202b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6215e76c431SBarry Smith d = b*b - 3*a*initslope; 6225e76c431SBarry Smith if (d < 0.0) d = 0.0; 6235e76c431SBarry Smith if (a == 0.0) { 6245e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6255e76c431SBarry Smith } else { 6265e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6275e76c431SBarry Smith } 6285e76c431SBarry Smith lambdaprev = lambda; 6295e76c431SBarry Smith gnormprev = *gnorm; 630329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 631329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6325e76c431SBarry Smith else lambda = lambdatemp; 633aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 634275d25c3SBarry Smith clambda = lambda; 635e68848bdSBarry Smith ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 6365e42470aSBarry Smith #else 637e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 6385e42470aSBarry Smith #endif 63943e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 640ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 641ae15b995SBarry 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); 642ed98166eSMatthew Knepley *flag = PETSC_FALSE; 64343e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 644ed98166eSMatthew Knepley break; 645ed98166eSMatthew Knepley } 64619717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 647d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 64819717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 64919717074SBarry Smith } 650d5e45103SBarry Smith CHKERRQ(ierr); 651cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 65243e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6535ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 654ae15b995SBarry Smith ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 65543e71028SBarry Smith break; 6562b022350SLois Curfman McInnes } else { 657ae15b995SBarry Smith ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 6585e76c431SBarry Smith } 6595e76c431SBarry Smith count++; 6605e76c431SBarry Smith } 6618229bfc2SMatthew Knepley if (count >= 10000) { 662cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6638229bfc2SMatthew Knepley } 664ad922ac9SBarry Smith theend1: 6658f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6663c632250SBarry Smith if (neP->postcheckstep && *flag) { 6673c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 6683c632250SBarry Smith if (changed_y) { 66979f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 6703c632250SBarry Smith } 6713c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 6723c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 673d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 67419717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 67519717074SBarry Smith } 676d5e45103SBarry Smith CHKERRQ(ierr); 6778f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 67843e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6793c632250SBarry Smith ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 6808f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6813c632250SBarry Smith ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 6828f99978dSLois Curfman McInnes } 6838f99978dSLois Curfman McInnes } 684d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6853a40ed3dSBarry Smith PetscFunctionReturn(0); 6865e76c431SBarry Smith } 68704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6884a2ae208SSatish Balay #undef __FUNCT__ 68917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic" 6904b828684SBarry Smith /*@C 69117bae607SBarry Smith SNESLineSearchQuadratic - Performs a quadratic line search. 6925e76c431SBarry Smith 693c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 694c7afd0dbSLois Curfman McInnes 6955e42470aSBarry Smith Input Parameters: 696c7afd0dbSLois Curfman McInnes + snes - the SNES context 697af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6985e42470aSBarry Smith . x - current iterate 6995e42470aSBarry Smith . f - residual evaluated at x 7003c632250SBarry Smith . y - search direction 7015e42470aSBarry Smith . w - work vector 702c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7035e42470aSBarry Smith 704c4a48953SLois Curfman McInnes Output Parameters: 7057f3332b4SBarry Smith + g - residual evaluated at new iterate w 7067f3332b4SBarry Smith . w - new iterate (x + alpha*y) 7075e42470aSBarry Smith . gnorm - 2-norm of g 7085e42470aSBarry Smith . ynorm - 2-norm of search length 709a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 710fee21e36SBarry Smith 711ce9499c7SBarry Smith Options Database Keys: 712ce9499c7SBarry Smith + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 713ce9499c7SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 714ce9499c7SBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 715ce9499c7SBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 716ce9499c7SBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 717ce9499c7SBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 7185e42470aSBarry Smith Notes: 71917bae607SBarry Smith Use SNESLineSearchSet() to set this routine within the SNESLS method. 7205e42470aSBarry Smith 72136851e7fSLois Curfman McInnes Level: advanced 72236851e7fSLois Curfman McInnes 723f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 724f59ffdedSLois Curfman McInnes 72517bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 7265e42470aSBarry Smith @*/ 72717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 7285e76c431SBarry Smith { 729406556e6SLois Curfman McInnes /* 730406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 731406556e6SLois Curfman McInnes minimization problem: 732406556e6SLois Curfman McInnes min z(x): R^n -> R, 733406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 734406556e6SLois Curfman McInnes */ 735275d25c3SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength; 736aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 737ea709b57SSatish Balay PetscScalar cinitslope,clambda; 7386b5873e3SBarry Smith #endif 739dfbe8321SBarry Smith PetscErrorCode ierr; 740a7cc72afSBarry Smith PetscInt count; 7414b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 74279f36460SBarry Smith PetscScalar scale; 7433c632250SBarry Smith PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 7445e76c431SBarry Smith 7453a40ed3dSBarry Smith PetscFunctionBegin; 746d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 747a7cc72afSBarry Smith *flag = PETSC_TRUE; 7485e42470aSBarry Smith alpha = neP->alpha; 7495e42470aSBarry Smith maxstep = neP->maxstep; 7505e42470aSBarry Smith steptol = neP->steptol; 7515e76c431SBarry Smith 7523505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 75363b9fa5eSBarry Smith if (*ynorm == 0.0) { 754ae15b995SBarry Smith ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 755b37302c6SLois Curfman McInnes *gnorm = fnorm; 756e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 757b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 758a7cc72afSBarry Smith *flag = PETSC_FALSE; 759ad922ac9SBarry Smith goto theend2; 76094a9d846SBarry Smith } 7615e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 7625e42470aSBarry Smith scale = maxstep/(*ynorm); 7632dcb1b2aSMatthew Knepley ierr = VecScale(y,scale);CHKERRQ(ierr); 7645e42470aSBarry Smith *ynorm = maxstep; 7655e76c431SBarry Smith } 766ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 7675ca10a37SBarry Smith minlambda = steptol/rellength; 768a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 769aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 770a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 771329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7725e42470aSBarry Smith #else 773a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7745e42470aSBarry Smith #endif 7755e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7765e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7775e42470aSBarry Smith 778e68848bdSBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 77943e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 780ae15b995SBarry Smith ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 78143e71028SBarry Smith *flag = PETSC_FALSE; 78243e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 78343e71028SBarry Smith goto theend2; 78443e71028SBarry Smith } 78519717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 786d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 78719717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 78819717074SBarry Smith } 789d5e45103SBarry Smith CHKERRQ(ierr); 790cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 79143e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 7925d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 793ae15b995SBarry Smith ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 794ad922ac9SBarry Smith goto theend2; 7955e42470aSBarry Smith } 7965e42470aSBarry Smith 7975e42470aSBarry Smith /* Fit points with quadratic */ 798313b5538SBarry Smith lambda = 1.0; 7995e42470aSBarry Smith count = 1; 8005ca10a37SBarry Smith while (PETSC_TRUE) { 8015e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 802ae15b995SBarry Smith ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 803ae15b995SBarry Smith ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 804e68848bdSBarry Smith ierr = VecCopy(x,w);CHKERRQ(ierr); 80543e71028SBarry Smith *flag = PETSC_FALSE; 80643e71028SBarry Smith break; 8075e42470aSBarry Smith } 808a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 809329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 810329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 811329f5518SBarry Smith else lambda = lambdatemp; 812275d25c3SBarry Smith 813aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 814e68848bdSBarry Smith clambda = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 8155e42470aSBarry Smith #else 816e68848bdSBarry Smith ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 8175e42470aSBarry Smith #endif 81843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 819ae15b995SBarry Smith ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 820ae15b995SBarry 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); 821ed98166eSMatthew Knepley *flag = PETSC_FALSE; 82243e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 823ed98166eSMatthew Knepley break; 824ed98166eSMatthew Knepley } 82519717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 826d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 82719717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 82819717074SBarry Smith } 829d5e45103SBarry Smith CHKERRQ(ierr); 830cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 83143e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8325ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 833ae15b995SBarry Smith ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 83406259719SBarry Smith break; 8355e42470aSBarry Smith } 8365e42470aSBarry Smith count++; 8375e42470aSBarry Smith } 838ad922ac9SBarry Smith theend2: 8398f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8403c632250SBarry Smith if (neP->postcheckstep) { 8413c632250SBarry Smith ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 8423c632250SBarry Smith if (changed_y) { 84379f36460SBarry Smith ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 8443c632250SBarry Smith } 8453c632250SBarry Smith if (changed_y || changed_w) { /* recompute the function if the step has changed */ 8463c632250SBarry Smith ierr = SNESComputeFunction(snes,w,g); 847d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 84819717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 84919717074SBarry Smith } 850d5e45103SBarry Smith CHKERRQ(ierr); 8518f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 8523c632250SBarry Smith ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 8538f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 8543c632250SBarry Smith ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 8553c632250SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8568f99978dSLois Curfman McInnes } 8578f99978dSLois Curfman McInnes } 858d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8593a40ed3dSBarry Smith PetscFunctionReturn(0); 8605e76c431SBarry Smith } 8612343ba6eSBarry Smith 86204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8634a2ae208SSatish Balay #undef __FUNCT__ 86417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet" 865c9e6a524SLois Curfman McInnes /*@C 86617bae607SBarry Smith SNESLineSearchSet - Sets the line search routine to be used 8674b27c08aSLois Curfman McInnes by the method SNESLS. 8685e76c431SBarry Smith 8695e76c431SBarry Smith Input Parameters: 870c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 871af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 872c7afd0dbSLois Curfman McInnes - func - pointer to int function 8735e76c431SBarry Smith 874fee21e36SBarry Smith Collective on SNES 875fee21e36SBarry Smith 876c4a48953SLois Curfman McInnes Available Routines: 87717bae607SBarry Smith + SNESLineSearchCubic() - default line search 87817bae607SBarry Smith . SNESLineSearchQuadratic() - quadratic line search 87917bae607SBarry Smith . SNESLineSearchNo() - the full Newton step (actually not a line search) 88017bae607SBarry Smith - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8815e76c431SBarry Smith 882c4a48953SLois Curfman McInnes Options Database Keys: 8834b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8844b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 8854b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 8864b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 8873304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 8883304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 889c4a48953SLois Curfman McInnes 8905e76c431SBarry Smith Calling sequence of func: 891bd208895SLois Curfman McInnes .vb 892af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 893a7cc72afSBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 894bd208895SLois Curfman McInnes .ve 8955e76c431SBarry Smith 8965e76c431SBarry Smith Input parameters for func: 897c7afd0dbSLois Curfman McInnes + snes - nonlinear context 898af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8995e76c431SBarry Smith . x - current iterate 9005e76c431SBarry Smith . f - residual evaluated at x 9013c632250SBarry Smith . y - search direction 9025e76c431SBarry Smith . w - work vector 903c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 9045e76c431SBarry Smith 9055e76c431SBarry Smith Output parameters for func: 906c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 9073c632250SBarry Smith . w - new iterate 9085e76c431SBarry Smith . gnorm - 2-norm of g 9095e76c431SBarry Smith . ynorm - 2-norm of search length 910a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 911f59ffdedSLois Curfman McInnes 91236851e7fSLois Curfman McInnes Level: advanced 91336851e7fSLois Curfman McInnes 914f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 915f59ffdedSLois Curfman McInnes 91617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 91717bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 9185e76c431SBarry Smith @*/ 91917bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 9205e76c431SBarry Smith { 921a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 92282bf6240SBarry Smith 9233a40ed3dSBarry Smith PetscFunctionBegin; 92417bae607SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 92582bf6240SBarry Smith if (f) { 926af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 92782bf6240SBarry Smith } 9283a40ed3dSBarry Smith PetscFunctionReturn(0); 9295e76c431SBarry Smith } 9308e019c35SBarry Smith 931a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 93204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 933fb2e594dSBarry Smith EXTERN_C_BEGIN 9344a2ae208SSatish Balay #undef __FUNCT__ 93517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS" 93617bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 93782bf6240SBarry Smith { 93882bf6240SBarry Smith PetscFunctionBegin; 9394b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9404b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 94182bf6240SBarry Smith PetscFunctionReturn(0); 94282bf6240SBarry Smith } 943fb2e594dSBarry Smith EXTERN_C_END 94404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9454a2ae208SSatish Balay #undef __FUNCT__ 9463c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck" 947c8dd0c0dSLois Curfman McInnes /*@C 9483c632250SBarry Smith SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 9494b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 950c8dd0c0dSLois Curfman McInnes 951c8dd0c0dSLois Curfman McInnes Input Parameters: 952c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 9533c632250SBarry Smith . func - pointer to function 954c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 955c8dd0c0dSLois Curfman McInnes 956c8dd0c0dSLois Curfman McInnes Collective on SNES 957c8dd0c0dSLois Curfman McInnes 958c8dd0c0dSLois Curfman McInnes Calling sequence of func: 959c8dd0c0dSLois Curfman McInnes .vb 9603c632250SBarry Smith int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 961c8dd0c0dSLois Curfman McInnes .ve 962b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 963b1ae27eaSLois Curfman McInnes on failure. 964c8dd0c0dSLois Curfman McInnes 965c8dd0c0dSLois Curfman McInnes Input parameters for func: 966c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 967c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9683c632250SBarry Smith . x - previous iterate 9693c632250SBarry Smith . y - new search direction and length 9703c632250SBarry Smith - w - current candidate iterate 971c8dd0c0dSLois Curfman McInnes 972c8dd0c0dSLois Curfman McInnes Output parameters for func: 9733c632250SBarry Smith + y - search direction (possibly changed) 9743c632250SBarry Smith . w - current iterate (possibly modified) 9753c632250SBarry Smith . changed_y - indicates search direction was changed by this routine 9763c632250SBarry Smith - changed_w - indicates current iterate was changed by this routine 977c8dd0c0dSLois Curfman McInnes 978c8dd0c0dSLois Curfman McInnes Level: advanced 979c8dd0c0dSLois Curfman McInnes 9809e247f21SBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 9819e247f21SBarry Smith 9823c632250SBarry Smith Only one of changed_y and changed_w can be PETSC_TRUE 9833c632250SBarry Smith 9843c632250SBarry Smith On input w = x + y 9853c632250SBarry Smith 98617bae607SBarry Smith SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 987b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 988ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9898f99978dSLois Curfman McInnes 99017bae607SBarry Smith SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 991ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 992ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 993ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 9949e247f21SBarry Smith by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 995b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9968f99978dSLois Curfman McInnes 997c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 998c8dd0c0dSLois Curfman McInnes 99917bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1000c8dd0c0dSLois Curfman McInnes @*/ 10013c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 1002c8dd0c0dSLois Curfman McInnes { 10033c632250SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 1004c8dd0c0dSLois Curfman McInnes 1005c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10063c632250SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1007c8dd0c0dSLois Curfman McInnes if (f) { 1008c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1009c8dd0c0dSLois Curfman McInnes } 1010c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1011c8dd0c0dSLois Curfman McInnes } 10129c3ca13aSBarry Smith #undef __FUNCT__ 10139c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck" 10149c3ca13aSBarry Smith /*@C 10159c3ca13aSBarry Smith SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 10167e4bb74cSBarry Smith before the line search is called. 10179c3ca13aSBarry Smith 10189c3ca13aSBarry Smith Input Parameters: 10199c3ca13aSBarry Smith + snes - nonlinear context obtained from SNESCreate() 10209c3ca13aSBarry Smith . func - pointer to function 10219c3ca13aSBarry Smith - checkctx - optional user-defined context for use by step checking routine 10229c3ca13aSBarry Smith 10239c3ca13aSBarry Smith Collective on SNES 10249c3ca13aSBarry Smith 10259c3ca13aSBarry Smith Calling sequence of func: 10269c3ca13aSBarry Smith .vb 10271e3347e8SBarry Smith int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 10289c3ca13aSBarry Smith .ve 10299c3ca13aSBarry Smith where func returns an error code of 0 on success and a nonzero 10309c3ca13aSBarry Smith on failure. 10319c3ca13aSBarry Smith 10329c3ca13aSBarry Smith Input parameters for func: 10339c3ca13aSBarry Smith + snes - nonlinear context 10349c3ca13aSBarry Smith . checkctx - optional user-defined context for use by step checking routine 10359c3ca13aSBarry Smith . x - previous iterate 10369c3ca13aSBarry Smith - y - new search direction and length 10379c3ca13aSBarry Smith 10389c3ca13aSBarry Smith Output parameters for func: 10399c3ca13aSBarry Smith + y - search direction (possibly changed) 10409c3ca13aSBarry Smith - changed_y - indicates search direction was changed by this routine 10419c3ca13aSBarry Smith 10429c3ca13aSBarry Smith Level: advanced 10439c3ca13aSBarry Smith 10449c3ca13aSBarry Smith Notes: All line searches accept the new iterate computed by the line search checking routine. 10459c3ca13aSBarry Smith 10469c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine 10479c3ca13aSBarry Smith 10487e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 10499c3ca13aSBarry Smith @*/ 10509c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 10519c3ca13aSBarry Smith { 10529c3ca13aSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 10539c3ca13aSBarry Smith 10549c3ca13aSBarry Smith PetscFunctionBegin; 10559c3ca13aSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 10569c3ca13aSBarry Smith if (f) { 10579c3ca13aSBarry Smith ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 10589c3ca13aSBarry Smith } 10599c3ca13aSBarry Smith PetscFunctionReturn(0); 10609c3ca13aSBarry Smith } 10619c3ca13aSBarry Smith 1062c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 10639c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 10649c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1065c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 10664a2ae208SSatish Balay #undef __FUNCT__ 10673c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 10689c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1069c8dd0c0dSLois Curfman McInnes { 1070c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10713c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheckstep = func; 10723c632250SBarry Smith ((SNES_LS *)(snes->data))->postcheck = checkctx; 1073c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1074c8dd0c0dSLois Curfman McInnes } 1075c8dd0c0dSLois Curfman McInnes EXTERN_C_END 10769c3ca13aSBarry Smith 10779c3ca13aSBarry Smith EXTERN_C_BEGIN 10789c3ca13aSBarry Smith #undef __FUNCT__ 10799c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 10809c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 10819c3ca13aSBarry Smith { 10829c3ca13aSBarry Smith PetscFunctionBegin; 10839c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheckstep = func; 10849c3ca13aSBarry Smith ((SNES_LS *)(snes->data))->precheck = checkctx; 10859c3ca13aSBarry Smith PetscFunctionReturn(0); 10869c3ca13aSBarry Smith } 10879c3ca13aSBarry Smith EXTERN_C_END 1088329e7e40SMatthew Knepley 1089329e7e40SMatthew Knepley /* 10904b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 109104d965bbSLois Curfman McInnes 109204d965bbSLois Curfman McInnes Input Parameters: 109304d965bbSLois Curfman McInnes . SNES - the SNES context 109404d965bbSLois Curfman McInnes . viewer - visualization context 109504d965bbSLois Curfman McInnes 109604d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 109704d965bbSLois Curfman McInnes */ 10984a2ae208SSatish Balay #undef __FUNCT__ 10994b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 11006849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1101a935fc98SLois Curfman McInnes { 11024b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 11032fc52814SBarry Smith const char *cstr; 1104dfbe8321SBarry Smith PetscErrorCode ierr; 110532077d6dSBarry Smith PetscTruth iascii; 1106a935fc98SLois Curfman McInnes 11073a40ed3dSBarry Smith PetscFunctionBegin; 110832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 110932077d6dSBarry Smith if (iascii) { 111017bae607SBarry Smith if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 111117bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 111217bae607SBarry Smith else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 111319bcc07fSBarry Smith else cstr = "unknown"; 1114b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1115a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 11165cd90555SBarry Smith } else { 111779a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 111819bcc07fSBarry Smith } 11193a40ed3dSBarry Smith PetscFunctionReturn(0); 1120a935fc98SLois Curfman McInnes } 112104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 112204d965bbSLois Curfman McInnes /* 11234b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 112404d965bbSLois Curfman McInnes 112504d965bbSLois Curfman McInnes Input Parameter: 112604d965bbSLois Curfman McInnes . snes - the SNES context 112704d965bbSLois Curfman McInnes 112804d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 112904d965bbSLois Curfman McInnes */ 11304a2ae208SSatish Balay #undef __FUNCT__ 11314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 11326849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 11335e42470aSBarry Smith { 11344b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1135e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1136dfbe8321SBarry Smith PetscErrorCode ierr; 1137a7cc72afSBarry Smith PetscInt indx; 1138f1af5d2fSBarry Smith PetscTruth flg; 11395e42470aSBarry Smith 11403a40ed3dSBarry Smith PetscFunctionBegin; 1141b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 11424b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 11434b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 11444b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1145186905e3SBarry Smith 114617bae607SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 114725cdf11fSBarry Smith if (flg) { 114822e36657SBarry Smith switch (indx) { 1149b49fd9e1SBarry Smith case 0: 115017bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1151b49fd9e1SBarry Smith break; 1152b49fd9e1SBarry Smith case 1: 115317bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1154b49fd9e1SBarry Smith break; 1155b49fd9e1SBarry Smith case 2: 115617bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1157b49fd9e1SBarry Smith break; 1158b49fd9e1SBarry Smith case 3: 115917bae607SBarry Smith ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1160b49fd9e1SBarry Smith break; 11615e42470aSBarry Smith } 11625e42470aSBarry Smith } 1163b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11643a40ed3dSBarry Smith PetscFunctionReturn(0); 11655e42470aSBarry Smith } 116604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 11673f149594SLisandro Dalcin #undef __FUNCT__ 11683f149594SLisandro Dalcin #define __FUNCT__ "SNESConverged_LS" 11693f149594SLisandro Dalcin /*@C 11703f149594SLisandro Dalcin SNESConverged_LS - Monitors the convergence of the line search 11713f149594SLisandro Dalcin method SNESLS for solving systems of nonlinear equations (default). 11723f149594SLisandro Dalcin 11733f149594SLisandro Dalcin Collective on SNES 11743f149594SLisandro Dalcin 11753f149594SLisandro Dalcin Input Parameters: 11763f149594SLisandro Dalcin + snes - the SNES context 11773f149594SLisandro Dalcin . it - the iteration (0 indicates before any Newton steps) 11783f149594SLisandro Dalcin . xnorm - 2-norm of current iterate 11793f149594SLisandro Dalcin . pnorm - 2-norm of current step 11803f149594SLisandro Dalcin . fnorm - 2-norm of function at current iterate 11813f149594SLisandro Dalcin - dummy - unused context 11823f149594SLisandro Dalcin 11833f149594SLisandro Dalcin Output Parameter: 11843f149594SLisandro Dalcin . reason - one of 11853f149594SLisandro Dalcin $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 11863f149594SLisandro Dalcin $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 11873f149594SLisandro Dalcin $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 11883f149594SLisandro Dalcin $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 11893f149594SLisandro Dalcin $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 11903f149594SLisandro Dalcin $ SNES_CONVERGED_ITERATING - (otherwise), 11913f149594SLisandro Dalcin 11923f149594SLisandro Dalcin where 11933f149594SLisandro Dalcin + maxf - maximum number of function evaluations, 11943f149594SLisandro Dalcin set with SNESSetTolerances() 11953f149594SLisandro Dalcin . nfct - number of function evaluations, 11963f149594SLisandro Dalcin . abstol - absolute function norm tolerance, 11973f149594SLisandro Dalcin set with SNESSetTolerances() 11983f149594SLisandro Dalcin - rtol - relative function norm tolerance, set with SNESSetTolerances() 11993f149594SLisandro Dalcin 12003f149594SLisandro Dalcin Level: intermediate 12013f149594SLisandro Dalcin 12023f149594SLisandro Dalcin .keywords: SNES, nonlinear, default, converged, convergence 12033f149594SLisandro Dalcin 12043f149594SLisandro Dalcin .seealso: SNESSetConvergenceTest() 12053f149594SLisandro Dalcin @*/ 12063f149594SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_LS(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 12073f149594SLisandro Dalcin { 12083f149594SLisandro Dalcin PetscErrorCode ierr; 12093f149594SLisandro Dalcin PetscFunctionBegin; 12103f149594SLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 12113f149594SLisandro Dalcin PetscValidType(snes,1); 12123f149594SLisandro Dalcin PetscValidPointer(reason,6); 12133f149594SLisandro Dalcin ierr = SNESDefaultConverged(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 12143f149594SLisandro Dalcin PetscFunctionReturn(0); 12153f149594SLisandro Dalcin } 12163f149594SLisandro Dalcin /* -------------------------------------------------------------------------- */ 1217ebe3b25bSBarry Smith /*MC 1218ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 121904d965bbSLois Curfman McInnes 1220ebe3b25bSBarry Smith Options Database: 1221ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1222ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1223ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1224ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1225ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1226ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 122704d965bbSLois Curfman McInnes 1228ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 122904d965bbSLois Curfman McInnes 1230ee3001cbSBarry Smith Level: beginner 1231ee3001cbSBarry Smith 123217bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 123317bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 123417bae607SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1235ebe3b25bSBarry Smith 1236ebe3b25bSBarry Smith M*/ 1237fb2e594dSBarry Smith EXTERN_C_BEGIN 12384a2ae208SSatish Balay #undef __FUNCT__ 12394b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 124063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 12415e42470aSBarry Smith { 1242dfbe8321SBarry Smith PetscErrorCode ierr; 12434b27c08aSLois Curfman McInnes SNES_LS *neP; 12445e42470aSBarry Smith 12453a40ed3dSBarry Smith PetscFunctionBegin; 1246e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 1247e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 1248e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 1249e7788613SBarry Smith snes->ops->converged = SNESConverged_LS; 1250e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 1251e7788613SBarry Smith snes->ops->view = SNESView_LS; 12525baf8537SBarry Smith snes->nwork = 0; 12535e42470aSBarry Smith 125438f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 12555e42470aSBarry Smith snes->data = (void*)neP; 12565e42470aSBarry Smith neP->alpha = 1.e-4; 12575e42470aSBarry Smith neP->maxstep = 1.e8; 12585e42470aSBarry Smith neP->steptol = 1.e-12; 125917bae607SBarry Smith neP->LineSearch = SNESLineSearchCubic; 1260c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 12613c632250SBarry Smith neP->postcheckstep = PETSC_NULL; 12623c632250SBarry Smith neP->postcheck = PETSC_NULL; 12633c632250SBarry Smith neP->precheckstep = PETSC_NULL; 12643c632250SBarry Smith neP->precheck = PETSC_NULL; 126582bf6240SBarry Smith 1266557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1267557d3f75SLisandro Dalcin "SNESLineSearchSet_LS", 1268557d3f75SLisandro Dalcin SNESLineSearchSet_LS);CHKERRQ(ierr); 1269557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1270557d3f75SLisandro Dalcin "SNESLineSearchSetPostCheck_LS", 1271557d3f75SLisandro Dalcin SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1272557d3f75SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1273557d3f75SLisandro Dalcin "SNESLineSearchSetPreCheck_LS", 1274557d3f75SLisandro Dalcin SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 127582bf6240SBarry Smith 12763a40ed3dSBarry Smith PetscFunctionReturn(0); 12775e42470aSBarry Smith } 1278fb2e594dSBarry Smith EXTERN_C_END 127982bf6240SBarry Smith 128082bf6240SBarry Smith 128182bf6240SBarry Smith 1282