1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> 3d4c6564cSPatrick Farrell 4d4c6564cSPatrick Farrell typedef struct { 5d4c6564cSPatrick Farrell PetscReal norm_delta_x_prev; /* norm of previous update */ 6d4c6564cSPatrick Farrell PetscReal norm_bar_delta_x_prev; /* norm of previous bar update */ 7d4c6564cSPatrick Farrell PetscReal mu_curr; /* current local Lipschitz estimate */ 8d4c6564cSPatrick Farrell PetscReal lambda_prev; /* previous step length: for some reason SNESLineSearchGetLambda returns 1 instead of the previous step length */ 9d4c6564cSPatrick Farrell } SNESLineSearch_NLEQERR; 10d4c6564cSPatrick Farrell 1157a6bf86SPatrick Farrell static PetscBool NLEQERR_cited = PETSC_FALSE; 1257a6bf86SPatrick Farrell static const char NLEQERR_citation[] = "@book{deuflhard2011,\n" 13d4c6564cSPatrick Farrell " title = {Newton Methods for Nonlinear Problems},\n" 14d4c6564cSPatrick Farrell " author = {Peter Deuflhard},\n" 15d4c6564cSPatrick Farrell " volume = 35,\n" 16d4c6564cSPatrick Farrell " year = 2011,\n" 17d4c6564cSPatrick Farrell " isbn = {978-3-642-23898-7},\n" 18d4c6564cSPatrick Farrell " doi = {10.1007/978-3-642-23899-4},\n" 19d4c6564cSPatrick Farrell " publisher = {Springer-Verlag},\n" 20d4c6564cSPatrick Farrell " address = {Berlin, Heidelberg}\n}\n"; 21d4c6564cSPatrick Farrell 229371c9d4SSatish Balay static PetscErrorCode SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch) { 2370d8d27cSBarry Smith SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data; 247cbffc34SPatrick Farrell 2570d8d27cSBarry Smith PetscFunctionBegin; 267cbffc34SPatrick Farrell nleqerr->mu_curr = 0.0; 277cbffc34SPatrick Farrell nleqerr->norm_delta_x_prev = -1.0; 287cbffc34SPatrick Farrell nleqerr->norm_bar_delta_x_prev = -1.0; 297cbffc34SPatrick Farrell PetscFunctionReturn(0); 307cbffc34SPatrick Farrell } 317cbffc34SPatrick Farrell 329371c9d4SSatish Balay static PetscErrorCode SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch) { 33d4c6564cSPatrick Farrell PetscBool changed_y, changed_w; 34d4c6564cSPatrick Farrell Vec X, F, Y, W, G; 35d4c6564cSPatrick Farrell SNES snes; 36d4c6564cSPatrick Farrell PetscReal fnorm, xnorm, ynorm, gnorm, wnorm; 37d4c6564cSPatrick Farrell PetscReal lambda, minlambda, stol; 38d4c6564cSPatrick Farrell PetscViewer monitor; 397cbffc34SPatrick Farrell PetscInt max_its, count, snes_iteration; 40d4c6564cSPatrick Farrell PetscReal theta, mudash, lambdadash; 4170d8d27cSBarry Smith SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data; 42d4c6564cSPatrick Farrell KSPConvergedReason kspreason; 43d4c6564cSPatrick Farrell 44d4c6564cSPatrick Farrell PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(NLEQERR_citation, &NLEQERR_cited)); 46d4c6564cSPatrick Farrell 479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G)); 489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 499566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 509566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 519566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 529566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, NULL, NULL, NULL, NULL, &max_its)); 539566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL)); 54d4c6564cSPatrick Farrell 557cbffc34SPatrick Farrell /* reset the state of the Lipschitz estimates */ 569566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &snes_iteration)); 57*48a46eb9SPierre Jolivet if (!snes_iteration) PetscCall(SNESLineSearchReset_NLEQERR(linesearch)); 587cbffc34SPatrick Farrell 59d4c6564cSPatrick Farrell /* precheck */ 609566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y)); 619566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 62d4c6564cSPatrick Farrell 639566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &ynorm)); 649566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 659566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &ynorm)); 669566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 67d4c6564cSPatrick Farrell 6870d8d27cSBarry Smith /* Note: Y is *minus* the Newton step. For whatever reason PETSc doesn't solve with the minus on the RHS. */ 69d4c6564cSPatrick Farrell 70d4c6564cSPatrick Farrell if (ynorm == 0.0) { 71d4c6564cSPatrick Farrell if (monitor) { 729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Initial direction and size is 0\n")); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 75d4c6564cSPatrick Farrell } 769566063dSJacob Faibussowitsch PetscCall(VecCopy(X, W)); 779566063dSJacob Faibussowitsch PetscCall(VecCopy(F, G)); 789566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm)); 799566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 80d4c6564cSPatrick Farrell PetscFunctionReturn(0); 81d4c6564cSPatrick Farrell } 82d4c6564cSPatrick Farrell 83d4c6564cSPatrick Farrell /* At this point, we've solved the Newton system for delta_x, and we assume that 84d4c6564cSPatrick Farrell its norm is greater than the solution tolerance (otherwise we wouldn't be in 85d4c6564cSPatrick Farrell here). So let's go ahead and estimate the Lipschitz constant. 86d4c6564cSPatrick Farrell 87d4c6564cSPatrick Farrell W contains bar_delta_x_prev at this point. */ 88d4c6564cSPatrick Farrell 89d4c6564cSPatrick Farrell if (monitor) { 909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: norm of Newton step: %14.12e\n", (double)ynorm)); 929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 93d4c6564cSPatrick Farrell } 94d4c6564cSPatrick Farrell 95d4c6564cSPatrick Farrell /* this needs information from a previous iteration, so can't do it on the first one */ 96d4c6564cSPatrick Farrell if (nleqerr->norm_delta_x_prev > 0 && nleqerr->norm_bar_delta_x_prev > 0) { 979566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, +1.0, Y, W)); /* bar_delta_x - delta_x; +1 because Y is -delta_x */ 989566063dSJacob Faibussowitsch PetscCall(VecNormBegin(G, NORM_2, &gnorm)); 999566063dSJacob Faibussowitsch PetscCall(VecNormEnd(G, NORM_2, &gnorm)); 100d4c6564cSPatrick Farrell 101d4c6564cSPatrick Farrell nleqerr->mu_curr = nleqerr->lambda_prev * (nleqerr->norm_delta_x_prev * nleqerr->norm_bar_delta_x_prev) / (gnorm * ynorm); 102d4c6564cSPatrick Farrell lambda = PetscMin(1.0, nleqerr->mu_curr); 103d4c6564cSPatrick Farrell 104d4c6564cSPatrick Farrell if (monitor) { 1059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: Lipschitz estimate: %14.12e; lambda: %14.12e\n", (double)nleqerr->mu_curr, (double)lambda)); 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 108d4c6564cSPatrick Farrell } 109bf8ddb9bSPatrick Farrell } else { 110d4c6564cSPatrick Farrell lambda = linesearch->damping; 111d4c6564cSPatrick Farrell } 112d4c6564cSPatrick Farrell 113d4c6564cSPatrick Farrell /* The main while loop of the algorithm. 114d4c6564cSPatrick Farrell At the end of this while loop, G should have the accepted new X in it. */ 115d4c6564cSPatrick Farrell 116d4c6564cSPatrick Farrell count = 0; 11759ea5d13SPatrick Farrell while (PETSC_TRUE) { 118d4c6564cSPatrick Farrell if (monitor) { 1199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 12063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: entering iteration with lambda: %14.12e\n", (double)lambda)); 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 122d4c6564cSPatrick Farrell } 123d4c6564cSPatrick Farrell 124d4c6564cSPatrick Farrell /* Check that we haven't performed too many iterations */ 125d4c6564cSPatrick Farrell count += 1; 126bf8ddb9bSPatrick Farrell if (count >= max_its) { 127d4c6564cSPatrick Farrell if (monitor) { 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: maximum iterations reached\n")); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 131d4c6564cSPatrick Farrell } 1329566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 133d4c6564cSPatrick Farrell PetscFunctionReturn(0); 134d4c6564cSPatrick Farrell } 135d4c6564cSPatrick Farrell 136d4c6564cSPatrick Farrell /* Now comes the Regularity Test. */ 137d4c6564cSPatrick Farrell if (lambda <= minlambda) { 138d4c6564cSPatrick Farrell /* This isn't what is suggested by Deuflhard, but it works better in my experience */ 139d4c6564cSPatrick Farrell if (monitor) { 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambda has reached lambdamin, taking full Newton step\n")); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 143d4c6564cSPatrick Farrell } 144d4c6564cSPatrick Farrell lambda = 1.0; 1459566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -lambda, Y, X)); 1467cbffc34SPatrick Farrell 1477cbffc34SPatrick Farrell /* and clean up the state for next time */ 1489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset_NLEQERR(linesearch)); 14970d8d27cSBarry Smith /* 15070d8d27cSBarry Smith The clang static analyzer detected a problem here; once the loop is broken the values 15170d8d27cSBarry Smith nleqerr->norm_delta_x_prev = ynorm; 15270d8d27cSBarry Smith nleqerr->norm_bar_delta_x_prev = wnorm; 15370d8d27cSBarry Smith are set, but wnorm has not even been computed. 15470d8d27cSBarry Smith I don't know if this is the correct fix but by setting ynorm and wnorm to -1.0 at 15570d8d27cSBarry Smith least the linesearch object is kept in the state set by the SNESLineSearchReset_NLEQERR() call above 15670d8d27cSBarry Smith */ 15770d8d27cSBarry Smith ynorm = wnorm = -1.0; 158d4c6564cSPatrick Farrell break; 159d4c6564cSPatrick Farrell } 160d4c6564cSPatrick Farrell 161d4c6564cSPatrick Farrell /* Compute new trial iterate */ 1629566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -lambda, Y, X)); 1639566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); 164d4c6564cSPatrick Farrell 165d4c6564cSPatrick Farrell /* Solve linear system for bar_delta_x_curr: old Jacobian, new RHS. Note absence of minus sign, compared to Deuflhard, in keeping with PETSc convention */ 1669566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, G, W)); 1679566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 168*48a46eb9SPierre Jolivet if (kspreason < 0) PetscCall(PetscInfo(snes, "Solution for \\bar{delta x}^{k+1} failed.")); 169d4c6564cSPatrick Farrell 170d4c6564cSPatrick Farrell /* W now contains -bar_delta_x_curr. */ 171d4c6564cSPatrick Farrell 1729566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 173d4c6564cSPatrick Farrell if (monitor) { 1749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: norm of simplified Newton update: %14.12e\n", (double)wnorm)); 1769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 177d4c6564cSPatrick Farrell } 178d4c6564cSPatrick Farrell 179d4c6564cSPatrick Farrell /* compute the monitoring quantities theta and mudash. */ 180d4c6564cSPatrick Farrell 181d4c6564cSPatrick Farrell theta = wnorm / ynorm; 182d4c6564cSPatrick Farrell 1839566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -(1.0 - lambda), Y, W)); 1849566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); 185d4c6564cSPatrick Farrell 186d4c6564cSPatrick Farrell mudash = (0.5 * ynorm * lambda * lambda) / gnorm; 187d4c6564cSPatrick Farrell 188d4c6564cSPatrick Farrell /* Check for termination of the linesearch */ 189d4c6564cSPatrick Farrell if (theta >= 1.0) { 190d4c6564cSPatrick Farrell /* need to go around again with smaller lambda */ 191d4c6564cSPatrick Farrell if (monitor) { 1929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel)); 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: monotonicity check failed, ratio: %14.12e\n", (double)theta)); 1949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel)); 195d4c6564cSPatrick Farrell } 196d4c6564cSPatrick Farrell lambda = PetscMin(mudash, 0.5 * lambda); 197d4c6564cSPatrick Farrell lambda = PetscMax(lambda, minlambda); 198d4c6564cSPatrick Farrell /* continue through the loop, i.e. go back to regularity test */ 199bf8ddb9bSPatrick Farrell } else { 200d4c6564cSPatrick Farrell /* linesearch terminated */ 201d4c6564cSPatrick Farrell lambdadash = PetscMin(1.0, mudash); 202d4c6564cSPatrick Farrell 203d4c6564cSPatrick Farrell if (lambdadash == 1.0 && lambda == 1.0 && wnorm <= stol) { 204d4c6564cSPatrick Farrell /* store the updated state, X - Y - W, in G: 205d4c6564cSPatrick Farrell I need to keep W for the next linesearch */ 206ef46b1a6SStefano Zampini PetscCall(VecWAXPY(G, -1.0, Y, X)); 2079566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, -1.0, W)); 208d4c6564cSPatrick Farrell break; 209d4c6564cSPatrick Farrell } 210d4c6564cSPatrick Farrell 211d4c6564cSPatrick Farrell /* Deuflhard suggests to add the following: 212d4c6564cSPatrick Farrell else if (lambdadash >= 4.0 * lambda) { 213d4c6564cSPatrick Farrell lambda = lambdadash; 214d4c6564cSPatrick Farrell } 215d4c6564cSPatrick Farrell to continue through the loop, i.e. go back to regularity test. 216d4c6564cSPatrick Farrell I deliberately exclude this, as I have practical experience of this 217d4c6564cSPatrick Farrell getting stuck in infinite loops (on e.g. an Allen--Cahn problem). */ 218d4c6564cSPatrick Farrell 219d4c6564cSPatrick Farrell else { 220d4c6564cSPatrick Farrell /* accept iterate without adding on, i.e. don't use bar_delta_x; 221d4c6564cSPatrick Farrell again, I need to keep W for the next linesearch */ 2229566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -lambda, Y, X)); 223d4c6564cSPatrick Farrell break; 224d4c6564cSPatrick Farrell } 225d4c6564cSPatrick Farrell } 226d4c6564cSPatrick Farrell } 227d4c6564cSPatrick Farrell 2281baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, G)); 229d4c6564cSPatrick Farrell 230d4c6564cSPatrick Farrell /* W currently contains -bar_delta_u. Scale it so that it contains bar_delta_u. */ 2319566063dSJacob Faibussowitsch PetscCall(VecScale(W, -1.0)); 232d4c6564cSPatrick Farrell 233d4c6564cSPatrick Farrell /* postcheck */ 2349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, G, &changed_y, &changed_w)); 235d4c6564cSPatrick Farrell if (changed_y || changed_w) { 2369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_USER)); 2379566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Changing the search direction here doesn't make sense.\n")); 238d4c6564cSPatrick Farrell PetscFunctionReturn(0); 239d4c6564cSPatrick Farrell } 240d4c6564cSPatrick Farrell 241d4c6564cSPatrick Farrell /* copy the solution and information from this iteration over */ 242d4c6564cSPatrick Farrell nleqerr->norm_delta_x_prev = ynorm; 243d4c6564cSPatrick Farrell nleqerr->norm_bar_delta_x_prev = wnorm; 244d4c6564cSPatrick Farrell nleqerr->lambda_prev = lambda; 245d4c6564cSPatrick Farrell 2469566063dSJacob Faibussowitsch PetscCall(VecCopy(G, X)); 2479566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2489566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); 2499566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 2509566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 2519566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, (ynorm < 0 ? PETSC_INFINITY : ynorm))); 252d4c6564cSPatrick Farrell PetscFunctionReturn(0); 253d4c6564cSPatrick Farrell } 254d4c6564cSPatrick Farrell 2559371c9d4SSatish Balay PetscErrorCode SNESLineSearchView_NLEQERR(SNESLineSearch linesearch, PetscViewer viewer) { 256d4c6564cSPatrick Farrell PetscBool iascii; 257d4c6564cSPatrick Farrell SNESLineSearch_NLEQERR *nleqerr; 258d4c6564cSPatrick Farrell 259d4c6564cSPatrick Farrell PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 261d4c6564cSPatrick Farrell nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data; 262d4c6564cSPatrick Farrell if (iascii) { 2639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " NLEQ-ERR affine-covariant linesearch")); 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " current local Lipschitz estimate omega=%e\n", (double)nleqerr->mu_curr)); 265d4c6564cSPatrick Farrell } 266d4c6564cSPatrick Farrell PetscFunctionReturn(0); 267d4c6564cSPatrick Farrell } 268d4c6564cSPatrick Farrell 2699371c9d4SSatish Balay static PetscErrorCode SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch) { 270d4c6564cSPatrick Farrell PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(PetscFree(linesearch->data)); 272d4c6564cSPatrick Farrell PetscFunctionReturn(0); 273d4c6564cSPatrick Farrell } 274d4c6564cSPatrick Farrell 275d4c6564cSPatrick Farrell /*MC 276d4c6564cSPatrick Farrell SNESLINESEARCHNLEQERR - Error-oriented affine-covariant globalised Newton algorithm of Deuflhard (2011). 277d4c6564cSPatrick Farrell 278d4c6564cSPatrick Farrell This linesearch is intended for Newton-type methods which are affine covariant. Affine covariance 279d4c6564cSPatrick Farrell means that Newton's method will give the same iterations for F(x) = 0 and AF(x) = 0 for a nonsingular 280d4c6564cSPatrick Farrell matrix A. This is a fundamental property; the philosophy of this linesearch is that globalisations 281d4c6564cSPatrick Farrell of Newton's method should carefully preserve it. 282d4c6564cSPatrick Farrell 283d4c6564cSPatrick Farrell For a discussion of the theory behind this algorithm, see 284d4c6564cSPatrick Farrell 285d4c6564cSPatrick Farrell @book{deuflhard2011, 286d4c6564cSPatrick Farrell title={Newton Methods for Nonlinear Problems}, 287d4c6564cSPatrick Farrell author={Deuflhard, P.}, 288d4c6564cSPatrick Farrell volume={35}, 289d4c6564cSPatrick Farrell year={2011}, 290d4c6564cSPatrick Farrell publisher={Springer-Verlag}, 291d4c6564cSPatrick Farrell address={Berlin, Heidelberg} 292d4c6564cSPatrick Farrell } 293d4c6564cSPatrick Farrell 294d4c6564cSPatrick Farrell Pseudocode is given on page 148. 295d4c6564cSPatrick Farrell 296d4c6564cSPatrick Farrell Options Database Keys: 297d4c6564cSPatrick Farrell + -snes_linesearch_damping<1.0> - initial step length 298d4c6564cSPatrick Farrell - -snes_linesearch_minlambda<1e-12> - minimum step length allowed 299d4c6564cSPatrick Farrell 300d4c6564cSPatrick Farrell Contributed by Patrick Farrell <patrick.farrell@maths.ox.ac.uk> 301d4c6564cSPatrick Farrell 302d4c6564cSPatrick Farrell Level: advanced 303d4c6564cSPatrick Farrell 304db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 305d4c6564cSPatrick Farrell M*/ 3069371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch) { 307d4c6564cSPatrick Farrell SNESLineSearch_NLEQERR *nleqerr; 308d4c6564cSPatrick Farrell 309d4c6564cSPatrick Farrell PetscFunctionBegin; 310d4c6564cSPatrick Farrell linesearch->ops->apply = SNESLineSearchApply_NLEQERR; 311d4c6564cSPatrick Farrell linesearch->ops->destroy = SNESLineSearchDestroy_NLEQERR; 312d4c6564cSPatrick Farrell linesearch->ops->setfromoptions = NULL; 313d4c6564cSPatrick Farrell linesearch->ops->reset = SNESLineSearchReset_NLEQERR; 314d4c6564cSPatrick Farrell linesearch->ops->view = SNESLineSearchView_NLEQERR; 315d4c6564cSPatrick Farrell linesearch->ops->setup = NULL; 316d4c6564cSPatrick Farrell 3179566063dSJacob Faibussowitsch PetscCall(PetscNewLog(linesearch, &nleqerr)); 318d4c6564cSPatrick Farrell 319d4c6564cSPatrick Farrell linesearch->data = (void *)nleqerr; 320d4c6564cSPatrick Farrell linesearch->max_its = 40; 3217cbffc34SPatrick Farrell SNESLineSearchReset_NLEQERR(linesearch); 322d4c6564cSPatrick Farrell PetscFunctionReturn(0); 323d4c6564cSPatrick Farrell } 324