1 #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2 #include <petsc/private/snesimpl.h> 3 4 typedef struct { 5 PetscReal alpha; /* sufficient decrease parameter */ 6 } SNESLineSearch_BT; 7 8 /*@ 9 SNESLineSearchBTSetAlpha - Sets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` `SNESLineSearch` variant. 10 11 Input Parameters: 12 + linesearch - linesearch context 13 - alpha - The descent parameter 14 15 Level: intermediate 16 17 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTGetAlpha()` 18 @*/ 19 PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha) 20 { 21 SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data; 22 PetscBool isbt; 23 24 PetscFunctionBegin; 25 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 26 PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, SNESLINESEARCHBT, &isbt)); 27 if (isbt) bt->alpha = alpha; 28 PetscFunctionReturn(PETSC_SUCCESS); 29 } 30 31 /*@ 32 SNESLineSearchBTGetAlpha - Gets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` variant that was set with `SNESLineSearchBTSetAlpha()` 33 34 Input Parameter: 35 . linesearch - linesearch context 36 37 Output Parameter: 38 . alpha - The descent parameter 39 40 Level: intermediate 41 42 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchGetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTSetAlpha()` 43 @*/ 44 PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha) 45 { 46 SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data; 47 PetscBool isbt; 48 49 PetscFunctionBegin; 50 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 51 PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, SNESLINESEARCHBT, &isbt)); 52 PetscCheck(isbt, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Not for type %s", ((PetscObject)linesearch)->type_name); 53 *alpha = bt->alpha; 54 PetscFunctionReturn(PETSC_SUCCESS); 55 } 56 57 static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) 58 { 59 SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data; 60 PetscBool changed_y, changed_w; 61 Vec X, F, Y, W, G; 62 SNES snes; 63 PetscReal fnorm, xnorm, ynorm, gnorm; 64 PetscReal lambda, lambdatemp, lambdaprev, minlambda, initslope, alpha, stol; 65 PetscReal t1, t2, a, b, d; 66 PetscReal f; 67 PetscReal g, gprev; 68 PetscViewer monitor; 69 PetscInt max_it, count; 70 Mat jac; 71 SNESObjectiveFn *objective; 72 const char *const ordStr[] = {"Linear", "Quadratic", "Cubic"}; 73 74 PetscFunctionBegin; 75 PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G)); 76 PetscCall(SNESLineSearchGetNorms(linesearch, NULL, &fnorm, NULL)); 77 PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 78 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 79 PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 80 PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, NULL, NULL, NULL, NULL, &max_it)); 81 PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL)); 82 PetscCall(SNESGetObjective(snes, &objective, NULL)); 83 alpha = bt->alpha; 84 85 PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL)); 86 PetscCheck(jac || objective, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 87 88 PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y)); 89 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 90 91 PetscCall(VecNormBegin(Y, NORM_2, &ynorm)); 92 PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 93 PetscCall(VecNormEnd(Y, NORM_2, &ynorm)); 94 PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 95 96 if (ynorm == 0.0) { 97 if (monitor) { 98 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 99 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Initial direction and size is 0\n")); 100 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 101 } 102 PetscCall(VecCopy(X, W)); 103 PetscCall(VecCopy(F, G)); 104 PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm)); 105 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 109 /* if the SNES has an objective set, use that instead of the function value */ 110 if (objective) { 111 PetscCall(SNESComputeObjective(snes, X, &f)); 112 SNESLineSearchCheckObjectiveDomainError(snes, f); 113 } else { 114 f = 0.5 * PetscSqr(fnorm); 115 } 116 117 /* compute the initial slope */ 118 if (objective) { 119 /* slope comes from the function (assumed to be the gradient of the objective) */ 120 PetscCall(VecDotRealPart(Y, F, &initslope)); 121 } else { 122 /* slope comes from the normal equations */ 123 PetscCall(MatMult(jac, Y, W)); 124 PetscCall(VecDotRealPart(F, W, &initslope)); 125 if (initslope > 0.0) initslope = -initslope; 126 if (initslope == 0.0) initslope = -1.0; 127 } 128 129 /* repeatedly cut lambda until the norm or objective function is not infinity or NaN or lambda is too small */ 130 while (PETSC_TRUE) { 131 PetscCall(VecWAXPY(W, -lambda, Y, X)); 132 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 133 if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 134 PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while checking full step length!\n")); 135 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 136 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION)); 137 PetscFunctionReturn(PETSC_SUCCESS); 138 } 139 140 if (objective) { 141 PetscCall(SNESComputeObjective(snes, W, &g)); 142 } else { 143 PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 144 if (linesearch->ops->vinorm) { 145 gnorm = fnorm; 146 PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm)); 147 } else { 148 PetscCall(VecNorm(G, NORM_2, &gnorm)); 149 } 150 g = 0.5 * PetscSqr(gnorm); 151 } 152 PetscCall(SNESLineSearchMonitor(linesearch)); 153 154 if (!PetscIsInfOrNanReal(g)) break; 155 if (monitor) { 156 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 157 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: objective function at lambdas = %g is infinity or NaN, cutting lambda\n", (double)lambda)); 158 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 159 } 160 if (lambda <= minlambda) SNESCheckFunctionDomainError(snes, g); 161 lambda *= .5; 162 } 163 164 if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm)); 165 if (g <= f + lambda * alpha * initslope) { /* Sufficient reduction or step tolerance convergence */ 166 if (monitor) { 167 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 168 if (!objective) { 169 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm)); 170 } else { 171 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Using full step: old obj %14.12e new obj %14.12e\n", (double)f, (double)g)); 172 } 173 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 174 } 175 } else { 176 if (stol * xnorm > ynorm) { 177 /* Since the full step didn't give sufficient decrease and the step is tiny, exit */ 178 PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm)); 179 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 180 if (monitor) { 181 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 182 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n", (double)ynorm, (double)(stol * xnorm))); 183 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 184 } 185 PetscFunctionReturn(PETSC_SUCCESS); 186 } 187 /* Here to avoid -Wmaybe-uninitiliazed warnings */ 188 lambdaprev = lambda; 189 gprev = g; 190 if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR) { 191 /* Fit points with quadratic */ 192 lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope)); 193 lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda); 194 195 PetscCall(VecWAXPY(W, -lambda, Y, X)); 196 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 197 if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 198 PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while attempting quadratic backtracking! %" PetscInt_FMT " \n", snes->nfuncs)); 199 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 200 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION)); 201 PetscFunctionReturn(PETSC_SUCCESS); 202 } 203 if (objective) { 204 PetscCall(SNESComputeObjective(snes, W, &g)); 205 SNESLineSearchCheckObjectiveDomainError(snes, g); 206 } else { 207 PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 208 if (linesearch->ops->vinorm) { 209 gnorm = fnorm; 210 PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm)); 211 } else { 212 PetscCall(VecNorm(G, NORM_2, &gnorm)); 213 } 214 SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm); 215 g = 0.5 * PetscSqr(gnorm); 216 } 217 if (PetscIsInfOrNanReal(g)) { 218 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF)); 219 PetscCall(PetscInfo(snes, "Aborted due to infinity or NaN in function evaluation\n")); 220 PetscFunctionReturn(PETSC_SUCCESS); 221 } 222 if (monitor) { 223 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 224 if (!objective) { 225 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: gnorm after quadratic fit %14.12e\n", (double)gnorm)); 226 } else { 227 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: obj after quadratic fit %14.12e\n", (double)g)); 228 } 229 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 230 } 231 } 232 if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR && g <= f + lambda * alpha * initslope) { /* sufficient reduction */ 233 if (monitor) { 234 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 235 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda)); 236 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 237 } 238 } else { 239 for (count = 0; count < max_it; count++) { 240 if (lambda <= minlambda) { 241 if (monitor) { 242 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 243 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: unable to find good step length! After %" PetscInt_FMT " tries \n", count)); 244 if (!objective) { 245 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope)); 246 } else { 247 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope)); 248 } 249 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 250 } 251 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 252 PetscFunctionReturn(PETSC_SUCCESS); 253 } 254 if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 255 /* Fit points with cubic */ 256 t1 = g - f - lambda * initslope; 257 t2 = gprev - f - lambdaprev * initslope; 258 a = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev); 259 b = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev); 260 d = b * b - 3 * a * initslope; 261 if (d < 0.0) d = 0.0; 262 if (a == 0.0) lambdatemp = -initslope / (2.0 * b); 263 else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a); 264 } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 265 lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope)); 266 } else if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { /* Just backtrack */ 267 lambdatemp = .5 * lambda; 268 } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "Line search order %" PetscInt_FMT " for type bt", linesearch->order); 269 lambdaprev = lambda; 270 gprev = g; 271 272 lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda); 273 PetscCall(VecWAXPY(W, -lambda, Y, X)); 274 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 275 if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 276 PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count)); 277 if (!objective) PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm, (double)gnorm, (double)ynorm, (double)lambda, (double)initslope)); 278 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION)); 279 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 280 PetscFunctionReturn(PETSC_SUCCESS); 281 } 282 if (objective) { 283 PetscCall(SNESComputeObjective(snes, W, &g)); 284 SNESLineSearchCheckObjectiveDomainError(snes, g); 285 } else { 286 PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 287 if (linesearch->ops->vinorm) { 288 gnorm = fnorm; 289 PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm)); 290 } else { 291 PetscCall(VecNorm(G, NORM_2, &gnorm)); 292 } 293 SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm); 294 g = 0.5 * PetscSqr(gnorm); 295 } 296 if (g <= f + lambda * alpha * initslope) { /* is reduction enough? */ 297 if (monitor) { 298 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 299 if (!objective) { 300 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: %s step, current gnorm %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)gnorm, (double)lambda)); 301 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 302 } else { 303 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: %s step, obj %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)g, (double)lambda)); 304 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 305 } 306 } 307 break; 308 } else if (monitor) { 309 PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 310 if (!objective) { 311 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: %s step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)gnorm, (double)lambda)); 312 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 313 } else { 314 PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: %s step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)g, (double)lambda)); 315 PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 316 } 317 } 318 } 319 } 320 } 321 322 /* postcheck */ 323 PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 324 PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w)); 325 if (changed_y) { 326 if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X)); 327 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 328 } 329 if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 330 PetscCall((*linesearch->ops->snesfunc)(snes, W, G)); 331 if (linesearch->ops->vinorm) { 332 gnorm = fnorm; 333 PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm)); 334 } else { 335 PetscCall(VecNorm(G, NORM_2, &gnorm)); 336 } 337 SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm); 338 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 339 } 340 341 /* copy the solution over */ 342 PetscCall(VecCopy(W, X)); 343 PetscCall(VecCopy(G, F)); 344 PetscCall(VecNorm(X, NORM_2, &xnorm)); 345 PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm)); 346 PetscFunctionReturn(PETSC_SUCCESS); 347 } 348 349 static PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 350 { 351 PetscBool isascii; 352 SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data; 353 354 PetscFunctionBegin; 355 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 356 if (isascii) { 357 if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 358 PetscCall(PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n")); 359 } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 360 PetscCall(PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n")); 361 } 362 PetscCall(PetscViewerASCIIPrintf(viewer, " alpha=%e\n", (double)bt->alpha)); 363 } 364 PetscFunctionReturn(PETSC_SUCCESS); 365 } 366 367 static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 368 { 369 PetscFunctionBegin; 370 PetscCall(PetscFree(linesearch->data)); 371 PetscFunctionReturn(PETSC_SUCCESS); 372 } 373 374 static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems PetscOptionsObject) 375 { 376 SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data; 377 378 PetscFunctionBegin; 379 PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options"); 380 PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL)); 381 PetscOptionsHeadEnd(); 382 PetscFunctionReturn(PETSC_SUCCESS); 383 } 384 385 /*MC 386 SNESLINESEARCHBT - Backtracking line search {cite}`dennis:83`. 387 388 This line search finds the minimum of a polynomial fitting either $1/2 ||F(x_k + \lambda Y_k)||_2^2$, 389 or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`. 390 If this fit does not satisfy the sufficient decrease conditions, the interval shrinks 391 and the fit is reattempted at most `max_it` times or until $\lambda$ is below `minlambda`. 392 393 Options Database Keys: 394 + -snes_linesearch_alpha <1e\-4> - slope descent parameter 395 . -snes_linesearch_damping <1.0> - initial `lambda` on entry to the line search 396 . -snes_linesearch_max_it <40> - maximum number of shrinking iterations in the line search 397 . -snes_linesearch_minlambda <1e\-12> - minimum `lambda` (scaling of solution update) allowed 398 - -snes_linesearch_order <3> - order of the polynomial fit, must be 1, 2, or 3. With order 1, it performs a simple backtracking without any curve fitting 399 400 Level: advanced 401 402 Note: 403 This line search will always produce a step that is less than or equal to, in length, the full step size. 404 405 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 406 M*/ 407 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 408 { 409 SNESLineSearch_BT *bt; 410 411 PetscFunctionBegin; 412 linesearch->ops->apply = SNESLineSearchApply_BT; 413 linesearch->ops->destroy = SNESLineSearchDestroy_BT; 414 linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 415 linesearch->ops->reset = NULL; 416 linesearch->ops->view = SNESLineSearchView_BT; 417 linesearch->ops->setup = NULL; 418 419 PetscCall(PetscNew(&bt)); 420 421 linesearch->data = (void *)bt; 422 linesearch->max_it = 40; 423 linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 424 bt->alpha = 1e-4; 425 PetscFunctionReturn(PETSC_SUCCESS); 426 } 427