1b9d635d7SJonas Heinzmann #include <petsc/private/linesearchimpl.h> 2b9d635d7SJonas Heinzmann #include <petsc/private/snesimpl.h> 3b9d635d7SJonas Heinzmann 4b9d635d7SJonas Heinzmann static PetscErrorCode SNESLineSearchApply_Bisection(SNESLineSearch linesearch) 5b9d635d7SJonas Heinzmann { 6b9d635d7SJonas Heinzmann PetscBool changed_y, changed_w; 7b9d635d7SJonas Heinzmann Vec X, F, Y, W, G; 8b9d635d7SJonas Heinzmann SNES snes; 9b9d635d7SJonas Heinzmann PetscReal ynorm; 10b9d635d7SJonas Heinzmann PetscReal lambda_left, lambda, lambda_right, lambda_old; 11*d5def619SJonas Heinzmann PetscScalar fty_left, fty, fty_initial; 12b9d635d7SJonas Heinzmann PetscViewer monitor; 13b9d635d7SJonas Heinzmann PetscReal rtol, atol, ltol; 14b9d635d7SJonas Heinzmann PetscInt it, max_its; 15b9d635d7SJonas Heinzmann 16b9d635d7SJonas Heinzmann PetscFunctionBegin; 17b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G)); 18b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 19b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 20b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetTolerances(linesearch, NULL, NULL, &rtol, &atol, <ol, &max_its)); 21b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 22b9d635d7SJonas Heinzmann 23b9d635d7SJonas Heinzmann /* pre-check */ 24b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y)); 25b9d635d7SJonas Heinzmann 26b9d635d7SJonas Heinzmann /* compute ynorm to normalize search direction */ 27b9d635d7SJonas Heinzmann PetscCall(VecNorm(Y, NORM_2, &ynorm)); 28b9d635d7SJonas Heinzmann 29b9d635d7SJonas Heinzmann /* initialize interval for bisection */ 30b9d635d7SJonas Heinzmann lambda_left = 0.0; 31b9d635d7SJonas Heinzmann lambda_right = lambda; 32b9d635d7SJonas Heinzmann 33b9d635d7SJonas Heinzmann /* compute fty at left end of interval */ 34*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 35*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, X, Y, &fty_left)); 36*d5def619SJonas Heinzmann } else { 37*d5def619SJonas Heinzmann PetscCall(VecDot(F, Y, &fty_left)); 38*d5def619SJonas Heinzmann } 39b9d635d7SJonas Heinzmann fty_initial = fty_left; 40b9d635d7SJonas Heinzmann 41b9d635d7SJonas Heinzmann /* compute fty at right end of interval (initial lambda) */ 42b9d635d7SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X)); 43b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 44b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 45b9d635d7SJonas Heinzmann if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 46b9d635d7SJonas Heinzmann PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations during line search!\n")); 47b9d635d7SJonas Heinzmann snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 48b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION)); 49b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 50b9d635d7SJonas Heinzmann } 51*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 52*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, G, W, Y, &fty)); 53*d5def619SJonas Heinzmann } else { 54*d5def619SJonas Heinzmann PetscCall(VecDot(G, Y, &fty)); 55*d5def619SJonas Heinzmann } 56b9d635d7SJonas Heinzmann /* check whether sign changes in interval */ 57*d5def619SJonas Heinzmann if (!PetscIsInfOrNanScalar(fty) && (PetscRealPart(fty_left * fty) > 0.0)) { 58b9d635d7SJonas Heinzmann /* no change of sign: accept full step */ 59b9d635d7SJonas Heinzmann if (monitor) { 60b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 61*d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: sign of fty does not change in step interval, accepting full step\n")); 62b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 63b9d635d7SJonas Heinzmann } 64b9d635d7SJonas Heinzmann } else { 65b9d635d7SJonas Heinzmann /* change of sign: iteratively bisect interval */ 66b9d635d7SJonas Heinzmann lambda_old = 0.0; 67b9d635d7SJonas Heinzmann it = 0; 68b9d635d7SJonas Heinzmann 69b9d635d7SJonas Heinzmann while (PETSC_TRUE) { 70b9d635d7SJonas Heinzmann /* check for NaN or Inf */ 71*d5def619SJonas Heinzmann if (PetscIsInfOrNanScalar(fty)) { 72b9d635d7SJonas Heinzmann if (monitor) { 73b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 74b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search fty is NaN or Inf!\n")); 75b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 76b9d635d7SJonas Heinzmann } 77b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF)); 78b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 79b9d635d7SJonas Heinzmann break; 80b9d635d7SJonas Heinzmann } 81b9d635d7SJonas Heinzmann 82b9d635d7SJonas Heinzmann /* check absolute tolerance */ 83*d5def619SJonas Heinzmann if (PetscAbsScalar(fty) <= atol * ynorm) { 84b9d635d7SJonas Heinzmann if (monitor) { 85b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 86*d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty)/||y|| = %g <= atol = %g\n", (double)(PetscAbsScalar(fty) / ynorm), (double)atol)); 87b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 88b9d635d7SJonas Heinzmann } 89b9d635d7SJonas Heinzmann break; 90b9d635d7SJonas Heinzmann } 91b9d635d7SJonas Heinzmann 92b9d635d7SJonas Heinzmann /* check relative tolerance */ 93*d5def619SJonas Heinzmann if (PetscAbsScalar(fty) / PetscAbsScalar(fty_initial) <= rtol) { 94b9d635d7SJonas Heinzmann if (monitor) { 95b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 96*d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty/fty_initial) = %g <= rtol = %g\n", (double)(PetscAbsScalar(fty) / PetscAbsScalar(fty_initial)), (double)rtol)); 97b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 98b9d635d7SJonas Heinzmann } 99b9d635d7SJonas Heinzmann break; 100b9d635d7SJonas Heinzmann } 101b9d635d7SJonas Heinzmann 102b9d635d7SJonas Heinzmann /* check maximum number of iterations */ 103b9d635d7SJonas Heinzmann if (it > max_its) { 104b9d635d7SJonas Heinzmann if (monitor) { 105b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 106b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: maximum iterations reached\n")); 107b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 108b9d635d7SJonas Heinzmann } 109b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 110b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 111b9d635d7SJonas Heinzmann break; 112b9d635d7SJonas Heinzmann } 113b9d635d7SJonas Heinzmann 114b9d635d7SJonas Heinzmann /* check change of lambda tolerance */ 115b9d635d7SJonas Heinzmann if (PetscAbsReal(lambda - lambda_old) < ltol) { 116b9d635d7SJonas Heinzmann if (monitor) { 117b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 118b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(dlambda) = %g < ltol = %g\n", (double)PetscAbsReal(lambda - lambda_old), (double)ltol)); 119b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 120b9d635d7SJonas Heinzmann } 121b9d635d7SJonas Heinzmann break; 122b9d635d7SJonas Heinzmann } 123b9d635d7SJonas Heinzmann 124b9d635d7SJonas Heinzmann /* determine direction of bisection (not necessary for 0th iteration) */ 125b9d635d7SJonas Heinzmann if (it > 0) { 126*d5def619SJonas Heinzmann if (PetscRealPart(fty * fty_left) <= 0.0) { 127b9d635d7SJonas Heinzmann lambda_right = lambda; 128b9d635d7SJonas Heinzmann } else { 129b9d635d7SJonas Heinzmann lambda_left = lambda; 130b9d635d7SJonas Heinzmann /* also update fty_left for direction check in next iteration */ 131b9d635d7SJonas Heinzmann fty_left = fty; 132b9d635d7SJonas Heinzmann } 133b9d635d7SJonas Heinzmann } 134b9d635d7SJonas Heinzmann 135b9d635d7SJonas Heinzmann /* bisect interval */ 136b9d635d7SJonas Heinzmann lambda_old = lambda; 137b9d635d7SJonas Heinzmann lambda = 0.5 * (lambda_left + lambda_right); 138b9d635d7SJonas Heinzmann 139b9d635d7SJonas Heinzmann /* compute fty at new lambda */ 140b9d635d7SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X)); 141b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 142b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 143b9d635d7SJonas Heinzmann if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 144b9d635d7SJonas Heinzmann PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations during line search!\n")); 145b9d635d7SJonas Heinzmann snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 146b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION)); 147b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 148b9d635d7SJonas Heinzmann } 149*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 150*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, G, W, Y, &fty)); 151*d5def619SJonas Heinzmann } else { 152*d5def619SJonas Heinzmann PetscCall(VecDot(G, Y, &fty)); 153*d5def619SJonas Heinzmann } 154b9d635d7SJonas Heinzmann 155b9d635d7SJonas Heinzmann /* print iteration information */ 156b9d635d7SJonas Heinzmann if (monitor) { 157b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 158*d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " %3" PetscInt_FMT " Line search: fty/||y|| = %g, lambda = %g\n", it, (double)(PetscRealPart(fty) / ynorm), (double)lambda)); 159b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 160b9d635d7SJonas Heinzmann } 161b9d635d7SJonas Heinzmann 162b9d635d7SJonas Heinzmann /* count up iteration */ 163b9d635d7SJonas Heinzmann it++; 164b9d635d7SJonas Heinzmann } 165b9d635d7SJonas Heinzmann } 166b9d635d7SJonas Heinzmann 167b9d635d7SJonas Heinzmann /* post-check */ 168b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 169b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w)); 170b9d635d7SJonas Heinzmann if (changed_y) { 171b9d635d7SJonas Heinzmann if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X)); 172b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 173b9d635d7SJonas Heinzmann } 174b9d635d7SJonas Heinzmann 175b9d635d7SJonas Heinzmann /* update solution*/ 176b9d635d7SJonas Heinzmann PetscCall(VecCopy(W, X)); 177b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, X, F)); 178b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchComputeNorms(linesearch)); 179b9d635d7SJonas Heinzmann 180b9d635d7SJonas Heinzmann /* finalization */ 181b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 182b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 183b9d635d7SJonas Heinzmann } 184b9d635d7SJonas Heinzmann 185b9d635d7SJonas Heinzmann /*MC 186b9d635d7SJonas Heinzmann SNESLINESEARCHBISECTION - Bisection line search. 187*d5def619SJonas Heinzmann Similar to the critical point line search, `SNESLINESEARCHCP`, the bisection line search assumes that there exists some $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$. 188b9d635d7SJonas Heinzmann This line search seeks to find the root of the directional derivative along the search direction $F^T Y$ through bisection. 189b9d635d7SJonas Heinzmann 190b9d635d7SJonas Heinzmann Options Database Keys: 191b9d635d7SJonas Heinzmann + -snes_linesearch_max_it <50> - maximum number of iterations for the line search 192b9d635d7SJonas Heinzmann . -snes_linesearch_damping <1.0> - initial trial step length on entry to the line search 193*d5def619SJonas Heinzmann . -snes_linesearch_rtol <1e\-8> - relative tolerance for the directional derivative 194*d5def619SJonas Heinzmann . -snes_linesearch_atol <1e\-6> - absolute tolerance for the directional derivative 195*d5def619SJonas Heinzmann - -snes_linesearch_ltol <1e\-6> - minimum absolute change in lambda allowed 196b9d635d7SJonas Heinzmann 197b9d635d7SJonas Heinzmann Level: intermediate 198b9d635d7SJonas Heinzmann 199b9d635d7SJonas Heinzmann Note: 200b9d635d7SJonas Heinzmann This method does NOT use the objective function if it is provided with `SNESSetObjective()`. 201b9d635d7SJonas Heinzmann This line search will always give a step size in the interval [0, damping]. 202b9d635d7SJonas Heinzmann 203b9d635d7SJonas Heinzmann .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLINESEARCHCP` 204b9d635d7SJonas Heinzmann M*/ 205b9d635d7SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Bisection(SNESLineSearch linesearch) 206b9d635d7SJonas Heinzmann { 207b9d635d7SJonas Heinzmann PetscFunctionBegin; 208b9d635d7SJonas Heinzmann linesearch->ops->apply = SNESLineSearchApply_Bisection; 209b9d635d7SJonas Heinzmann linesearch->ops->destroy = NULL; 210b9d635d7SJonas Heinzmann linesearch->ops->setfromoptions = NULL; 211b9d635d7SJonas Heinzmann linesearch->ops->reset = NULL; 212b9d635d7SJonas Heinzmann linesearch->ops->view = NULL; 213b9d635d7SJonas Heinzmann linesearch->ops->setup = NULL; 214b9d635d7SJonas Heinzmann 215b9d635d7SJonas Heinzmann /* set default option values */ 216b9d635d7SJonas Heinzmann linesearch->max_its = 50; 217b9d635d7SJonas Heinzmann linesearch->rtol = 1e-8; 218b9d635d7SJonas Heinzmann linesearch->atol = 1e-6; 219b9d635d7SJonas Heinzmann linesearch->ltol = 1e-6; 220b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS); 221b9d635d7SJonas Heinzmann } 222