1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> 26a388c36SPeter Brune #include <petscsnes.h> 3bbd5d0b3SPeter Brune 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch) 5d71ae5a4SJacob Faibussowitsch { 6422a814eSBarry Smith PetscBool changed_y, changed_w; 76a388c36SPeter Brune Vec X, Y, F, W; 86a388c36SPeter Brune SNES snes; 9516fe3c3SPeter Brune PetscReal xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep; 10bbd5d0b3SPeter Brune PetscReal lambda, lambda_old, lambda_update, delLambda; 11b8ac21a2SPeter Brune PetscScalar fty, fty_init, fty_old, fty_mid1, fty_mid2, s; 126a388c36SPeter Brune PetscInt i, max_its; 136a388c36SPeter Brune PetscViewer monitor; 14bbd5d0b3SPeter Brune 15bbd5d0b3SPeter Brune PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL)); 179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 189566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 199566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 209566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its)); 219566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 229566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 236a388c36SPeter Brune 24bbd5d0b3SPeter Brune /* precheck */ 259566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y)); 26bbd5d0b3SPeter Brune lambda_old = 0.0; 2762893cf2SPeter Brune 289566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_old)); 2996b7d5c2SJed Brown if (PetscAbsScalar(fty_old) < atol * ynorm) { 309e0a2ef9SBarry Smith if (monitor) { 319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 3263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search terminated at initial point because dot(F,Y) = %g < atol*||y|| = %g\n", (double)PetscAbsScalar(fty_old), (double)(atol * ynorm))); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 349e0a2ef9SBarry Smith } 359566063dSJacob Faibussowitsch PetscCall(SNESSetConvergedReason(linesearch->snes, SNES_CONVERGED_FNORM_ABS)); 363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 379e0a2ef9SBarry Smith } 389e0a2ef9SBarry Smith 3962893cf2SPeter Brune fty_init = fty_old; 40c0b2db79SJed Brown 41c0b2db79SJed Brown for (i = 0; i < max_its; i++) { 42bbd5d0b3SPeter Brune /* compute the norm at lambda */ 43ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -lambda, Y, X)); 441baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 459566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 469566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty)); 47bbd5d0b3SPeter Brune 48bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 49b8ac21a2SPeter Brune 50b8ac21a2SPeter Brune /* check for convergence */ 51b8ac21a2SPeter Brune if (PetscAbsReal(delLambda) < steptol * lambda) break; 52b8ac21a2SPeter Brune if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 5396b7d5c2SJed Brown if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break; 546a388c36SPeter Brune if (monitor) { 559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambdas = [%g, %g], ftys = [%g, %g]\n", (double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old))); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 58bbd5d0b3SPeter Brune } 59bbd5d0b3SPeter Brune 60bbd5d0b3SPeter Brune /* compute the search direction */ 61b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { 62b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda; 63b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 64ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X)); 651baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 669566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 679566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid1)); 68b8ac21a2SPeter Brune s = (3. * fty - 4. * fty_mid1 + fty_old) / delLambda; 69b8ac21a2SPeter Brune } else { 70ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X)); 711baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 729566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 739566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid1)); 74ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -(lambda + 0.5 * (lambda - lambda_old)), Y, X)); 751baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 769566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 779566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid2)); 78b8ac21a2SPeter Brune s = (2. * fty_mid2 + 3. * fty - 6. * fty_mid1 + fty_old) / (3. * delLambda); 79b8ac21a2SPeter Brune } 80b8ac21a2SPeter Brune /* if the solve is going in the wrong direction, fix it */ 81b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s; 82f2f03862SBarry Smith if (s == 0.0) break; 83b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s); 84b8ac21a2SPeter Brune 85b8ac21a2SPeter Brune /* switch directions if we stepped out of bounds */ 86f5af7f23SKarl Rupp if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s); 87b8ac21a2SPeter Brune 88620fc8aaSSatish Balay if (PetscIsInfOrNanReal(lambda_update)) break; 89f5af7f23SKarl Rupp if (lambda_update > maxstep) break; 90bbd5d0b3SPeter Brune 91bbd5d0b3SPeter Brune /* compute the new state of the line search */ 92bbd5d0b3SPeter Brune lambda_old = lambda; 93bbd5d0b3SPeter Brune lambda = lambda_update; 94bbd5d0b3SPeter Brune fty_old = fty; 95bbd5d0b3SPeter Brune } 96bbd5d0b3SPeter Brune /* construct the solution */ 97ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -lambda, Y, X)); 981baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 99bbd5d0b3SPeter Brune /* postcheck */ 1009566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w)); 101*a42aa960SStefano Zampini if (changed_y && !changed_w) { 1029566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -lambda, Y)); 1031baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, X)); 104bbd5d0b3SPeter Brune } else { 1059566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 106bbd5d0b3SPeter Brune } 1079566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, X, F)); 108bbd5d0b3SPeter Brune 1099566063dSJacob Faibussowitsch PetscCall(SNESLineSearchComputeNorms(linesearch)); 1109566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 111bbd5d0b3SPeter Brune 1129566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 1136a388c36SPeter Brune 1146a388c36SPeter Brune if (monitor) { 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm)); 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 118bbd5d0b3SPeter Brune } 11948a46eb9SPierre Jolivet if (lambda <= steptol) PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 121bbd5d0b3SPeter Brune } 122bbd5d0b3SPeter Brune 123954494b2SJed Brown /*MC 1241a4f838cSPeter Brune SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some 125f6dfbefdSBarry Smith artificial G(x) for which the `SNESFunction` F(x) = grad G(x). Therefore, this line search seeks 126db609ea7SPeter Brune to find roots of dot(F, Y) via a secant method. 127cd7522eaSPeter Brune 128cd7522eaSPeter Brune Options Database Keys: 129eb23ba39SBarry Smith + -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda 130eb23ba39SBarry Smith . -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value 131eb23ba39SBarry Smith . -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor, default is 1.0 132eb23ba39SBarry Smith - -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed. 133cd7522eaSPeter Brune 134cd7522eaSPeter Brune Notes: 135f6dfbefdSBarry Smith This method does NOT use the objective function if it is provided with `SNESSetObjective()`. 136eb23ba39SBarry Smith 137f6dfbefdSBarry Smith This method is the preferred line search for `SNESQN` and `SNESNCG`. 138954494b2SJed Brown 139954494b2SJed Brown Level: advanced 140954494b2SJed Brown 141f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 142954494b2SJed Brown M*/ 143d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 144d71ae5a4SJacob Faibussowitsch { 145bbd5d0b3SPeter Brune PetscFunctionBegin; 146f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP; 1470298fd71SBarry Smith linesearch->ops->destroy = NULL; 1480298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1490298fd71SBarry Smith linesearch->ops->reset = NULL; 1500298fd71SBarry Smith linesearch->ops->view = NULL; 1510298fd71SBarry Smith linesearch->ops->setup = NULL; 152b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_LINEAR; 153a491bab8SPeter Brune 154a491bab8SPeter Brune linesearch->max_its = 1; 1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 156bbd5d0b3SPeter Brune } 157