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 28*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 29*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, X, Y, &fty_old)); 30*d5def619SJonas Heinzmann } else { 319566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_old)); 32*d5def619SJonas Heinzmann } 3396b7d5c2SJed Brown if (PetscAbsScalar(fty_old) < atol * ynorm) { 349e0a2ef9SBarry Smith if (monitor) { 359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 3663a3b9bcSJacob 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))); 379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 389e0a2ef9SBarry Smith } 399566063dSJacob Faibussowitsch PetscCall(SNESSetConvergedReason(linesearch->snes, SNES_CONVERGED_FNORM_ABS)); 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 419e0a2ef9SBarry Smith } 429e0a2ef9SBarry Smith 4362893cf2SPeter Brune fty_init = fty_old; 44c0b2db79SJed Brown 45c0b2db79SJed Brown for (i = 0; i < max_its; i++) { 46bbd5d0b3SPeter Brune /* compute the norm at lambda */ 47ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -lambda, Y, X)); 481baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 499566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 50*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 51*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty)); 52*d5def619SJonas Heinzmann } else { 539566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty)); 54*d5def619SJonas Heinzmann } 55bbd5d0b3SPeter Brune 56bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 57b8ac21a2SPeter Brune 58b8ac21a2SPeter Brune /* check for convergence */ 59b8ac21a2SPeter Brune if (PetscAbsReal(delLambda) < steptol * lambda) break; 60b8ac21a2SPeter Brune if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 6196b7d5c2SJed Brown if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break; 626a388c36SPeter Brune if (monitor) { 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 649566063dSJacob 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))); 659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 66bbd5d0b3SPeter Brune } 67bbd5d0b3SPeter Brune 68bbd5d0b3SPeter Brune /* compute the search direction */ 69b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { 70b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda; 71b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 72ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X)); 731baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 749566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 75*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 76*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid1)); 77*d5def619SJonas Heinzmann } else { 789566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid1)); 79*d5def619SJonas Heinzmann } 80b8ac21a2SPeter Brune s = (3. * fty - 4. * fty_mid1 + fty_old) / delLambda; 81b8ac21a2SPeter Brune } else { 82ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X)); 831baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 849566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 85*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 86*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid1)); 87*d5def619SJonas Heinzmann } else { 889566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid1)); 89*d5def619SJonas Heinzmann } 90ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -(lambda + 0.5 * (lambda - lambda_old)), Y, X)); 911baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 929566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 93*d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) { 94*d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid2)); 95*d5def619SJonas Heinzmann } else { 969566063dSJacob Faibussowitsch PetscCall(VecDot(F, Y, &fty_mid2)); 97*d5def619SJonas Heinzmann } 98b8ac21a2SPeter Brune s = (2. * fty_mid2 + 3. * fty - 6. * fty_mid1 + fty_old) / (3. * delLambda); 99b8ac21a2SPeter Brune } 100b8ac21a2SPeter Brune /* if the solve is going in the wrong direction, fix it */ 101b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s; 102f2f03862SBarry Smith if (s == 0.0) break; 103b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s); 104b8ac21a2SPeter Brune 105b8ac21a2SPeter Brune /* switch directions if we stepped out of bounds */ 106f5af7f23SKarl Rupp if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s); 107b8ac21a2SPeter Brune 108620fc8aaSSatish Balay if (PetscIsInfOrNanReal(lambda_update)) break; 109f5af7f23SKarl Rupp if (lambda_update > maxstep) break; 110bbd5d0b3SPeter Brune 111bbd5d0b3SPeter Brune /* compute the new state of the line search */ 112bbd5d0b3SPeter Brune lambda_old = lambda; 113bbd5d0b3SPeter Brune lambda = lambda_update; 114bbd5d0b3SPeter Brune fty_old = fty; 115bbd5d0b3SPeter Brune } 116bbd5d0b3SPeter Brune /* construct the solution */ 117ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -lambda, Y, X)); 1181baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 119bbd5d0b3SPeter Brune /* postcheck */ 1206b095a85SStefano Zampini PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 1219566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w)); 1226b095a85SStefano Zampini if (changed_y) { 1236b095a85SStefano Zampini if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X)); 1246b095a85SStefano Zampini if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 125bbd5d0b3SPeter Brune } 1266b095a85SStefano Zampini PetscCall(VecCopy(W, X)); 1279566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, X, F)); 128bbd5d0b3SPeter Brune 1299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchComputeNorms(linesearch)); 1309566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 131bbd5d0b3SPeter Brune 1326a388c36SPeter Brune if (monitor) { 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 136bbd5d0b3SPeter Brune } 13748a46eb9SPierre Jolivet if (lambda <= steptol) PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 139bbd5d0b3SPeter Brune } 140bbd5d0b3SPeter Brune 141954494b2SJed Brown /*MC 1421a4f838cSPeter Brune SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some 143*d5def619SJonas Heinzmann artificial $G(x)$ for which the `SNESFunctionFn` $ F(x) = grad G(x)$. Therefore, this line search seeks 14489dfbbd5SBarry Smith to find roots of $ F^T Y$ via a secant method. 145cd7522eaSPeter Brune 146cd7522eaSPeter Brune Options Database Keys: 147eb23ba39SBarry Smith + -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda 148eb23ba39SBarry Smith . -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value 14989dfbbd5SBarry Smith . -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor on entry to the line search, default is 1.0 150eb23ba39SBarry Smith - -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed. 151cd7522eaSPeter Brune 152420bcc1bSBarry Smith Level: advanced 153420bcc1bSBarry Smith 154cd7522eaSPeter Brune Notes: 155f6dfbefdSBarry Smith This method does NOT use the objective function if it is provided with `SNESSetObjective()`. 156eb23ba39SBarry Smith 157f6dfbefdSBarry Smith This method is the preferred line search for `SNESQN` and `SNESNCG`. 158954494b2SJed Brown 159b9d635d7SJonas Heinzmann .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLINESEARCHBISECTION` 160954494b2SJed Brown M*/ 161d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 162d71ae5a4SJacob Faibussowitsch { 163bbd5d0b3SPeter Brune PetscFunctionBegin; 164f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP; 1650298fd71SBarry Smith linesearch->ops->destroy = NULL; 1660298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1670298fd71SBarry Smith linesearch->ops->reset = NULL; 1680298fd71SBarry Smith linesearch->ops->view = NULL; 1690298fd71SBarry Smith linesearch->ops->setup = NULL; 170b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_LINEAR; 171a491bab8SPeter Brune 172a491bab8SPeter Brune linesearch->max_its = 1; 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 174bbd5d0b3SPeter Brune } 175