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