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 227cbffc34SPatrick Farrell static PetscErrorCode SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch) 237cbffc34SPatrick Farrell { 2470d8d27cSBarry Smith SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR*)linesearch->data; 257cbffc34SPatrick Farrell 2670d8d27cSBarry Smith PetscFunctionBegin; 277cbffc34SPatrick Farrell nleqerr->mu_curr = 0.0; 287cbffc34SPatrick Farrell nleqerr->norm_delta_x_prev = -1.0; 297cbffc34SPatrick Farrell nleqerr->norm_bar_delta_x_prev = -1.0; 307cbffc34SPatrick Farrell PetscFunctionReturn(0); 317cbffc34SPatrick Farrell } 327cbffc34SPatrick Farrell 33d4c6564cSPatrick Farrell static PetscErrorCode SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch) 34d4c6564cSPatrick Farrell { 35d4c6564cSPatrick Farrell PetscBool changed_y,changed_w; 36d4c6564cSPatrick Farrell Vec X,F,Y,W,G; 37d4c6564cSPatrick Farrell SNES snes; 38d4c6564cSPatrick Farrell PetscReal fnorm, xnorm, ynorm, gnorm, wnorm; 39d4c6564cSPatrick Farrell PetscReal lambda, minlambda, stol; 40d4c6564cSPatrick Farrell PetscViewer monitor; 417cbffc34SPatrick Farrell PetscInt max_its, count, snes_iteration; 42d4c6564cSPatrick Farrell PetscReal theta, mudash, lambdadash; 4370d8d27cSBarry Smith SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR*)linesearch->data; 44d4c6564cSPatrick Farrell KSPConvergedReason kspreason; 45d4c6564cSPatrick Farrell 46d4c6564cSPatrick Farrell PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(NLEQERR_citation, &NLEQERR_cited)); 48d4c6564cSPatrick Farrell 499566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G)); 509566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 519566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 529566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 539566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 549566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetTolerances(linesearch,&minlambda,NULL,NULL,NULL,NULL,&max_its)); 559566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL)); 56d4c6564cSPatrick Farrell 577cbffc34SPatrick Farrell /* reset the state of the Lipschitz estimates */ 589566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &snes_iteration)); 5970d8d27cSBarry Smith if (!snes_iteration) { 609566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset_NLEQERR(linesearch)); 617cbffc34SPatrick Farrell } 627cbffc34SPatrick Farrell 63d4c6564cSPatrick Farrell /* precheck */ 649566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPreCheck(linesearch,X,Y,&changed_y)); 659566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 66d4c6564cSPatrick Farrell 679566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &ynorm)); 689566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 699566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &ynorm)); 709566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 71d4c6564cSPatrick Farrell 7270d8d27cSBarry Smith /* Note: Y is *minus* the Newton step. For whatever reason PETSc doesn't solve with the minus on the RHS. */ 73d4c6564cSPatrick Farrell 74d4c6564cSPatrick Farrell if (ynorm == 0.0) { 75d4c6564cSPatrick Farrell if (monitor) { 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n")); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 79d4c6564cSPatrick Farrell } 809566063dSJacob Faibussowitsch PetscCall(VecCopy(X,W)); 819566063dSJacob Faibussowitsch PetscCall(VecCopy(F,G)); 829566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm)); 839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 84d4c6564cSPatrick Farrell PetscFunctionReturn(0); 85d4c6564cSPatrick Farrell } 86d4c6564cSPatrick Farrell 87d4c6564cSPatrick Farrell /* At this point, we've solved the Newton system for delta_x, and we assume that 88d4c6564cSPatrick Farrell its norm is greater than the solution tolerance (otherwise we wouldn't be in 89d4c6564cSPatrick Farrell here). So let's go ahead and estimate the Lipschitz constant. 90d4c6564cSPatrick Farrell 91d4c6564cSPatrick Farrell W contains bar_delta_x_prev at this point. */ 92d4c6564cSPatrick Farrell 93d4c6564cSPatrick Farrell if (monitor) { 949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: norm of Newton step: %14.12e\n", (double) ynorm)); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 97d4c6564cSPatrick Farrell } 98d4c6564cSPatrick Farrell 99d4c6564cSPatrick Farrell /* this needs information from a previous iteration, so can't do it on the first one */ 100d4c6564cSPatrick Farrell if (nleqerr->norm_delta_x_prev > 0 && nleqerr->norm_bar_delta_x_prev > 0) { 1019566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, +1.0, Y, W)); /* bar_delta_x - delta_x; +1 because Y is -delta_x */ 1029566063dSJacob Faibussowitsch PetscCall(VecNormBegin(G, NORM_2, &gnorm)); 1039566063dSJacob Faibussowitsch PetscCall(VecNormEnd(G, NORM_2, &gnorm)); 104d4c6564cSPatrick Farrell 105d4c6564cSPatrick Farrell nleqerr->mu_curr = nleqerr->lambda_prev * (nleqerr->norm_delta_x_prev * nleqerr->norm_bar_delta_x_prev) / (gnorm * ynorm); 106d4c6564cSPatrick Farrell lambda = PetscMin(1.0, nleqerr->mu_curr); 107d4c6564cSPatrick Farrell 108d4c6564cSPatrick Farrell if (monitor) { 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: Lipschitz estimate: %14.12e; lambda: %14.12e\n", (double) nleqerr->mu_curr, (double) lambda)); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 112d4c6564cSPatrick Farrell } 113bf8ddb9bSPatrick Farrell } else { 114d4c6564cSPatrick Farrell lambda = linesearch->damping; 115d4c6564cSPatrick Farrell } 116d4c6564cSPatrick Farrell 117d4c6564cSPatrick Farrell /* The main while loop of the algorithm. 118d4c6564cSPatrick Farrell At the end of this while loop, G should have the accepted new X in it. */ 119d4c6564cSPatrick Farrell 120d4c6564cSPatrick Farrell count = 0; 12159ea5d13SPatrick Farrell while (PETSC_TRUE) { 122d4c6564cSPatrick Farrell if (monitor) { 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: entering iteration with lambda: %14.12e\n", lambda)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 126d4c6564cSPatrick Farrell } 127d4c6564cSPatrick Farrell 128d4c6564cSPatrick Farrell /* Check that we haven't performed too many iterations */ 129d4c6564cSPatrick Farrell count += 1; 130bf8ddb9bSPatrick Farrell if (count >= max_its) { 131d4c6564cSPatrick Farrell if (monitor) { 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: maximum iterations reached\n")); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 135d4c6564cSPatrick Farrell } 1369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT)); 137d4c6564cSPatrick Farrell PetscFunctionReturn(0); 138d4c6564cSPatrick Farrell } 139d4c6564cSPatrick Farrell 140d4c6564cSPatrick Farrell /* Now comes the Regularity Test. */ 141d4c6564cSPatrick Farrell if (lambda <= minlambda) { 142d4c6564cSPatrick Farrell /* This isn't what is suggested by Deuflhard, but it works better in my experience */ 143d4c6564cSPatrick Farrell if (monitor) { 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: lambda has reached lambdamin, taking full Newton step\n")); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 147d4c6564cSPatrick Farrell } 148d4c6564cSPatrick Farrell lambda = 1.0; 1499566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -lambda, Y, X)); 1507cbffc34SPatrick Farrell 1517cbffc34SPatrick Farrell /* and clean up the state for next time */ 1529566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset_NLEQERR(linesearch)); 15370d8d27cSBarry Smith /* 15470d8d27cSBarry Smith The clang static analyzer detected a problem here; once the loop is broken the values 15570d8d27cSBarry Smith nleqerr->norm_delta_x_prev = ynorm; 15670d8d27cSBarry Smith nleqerr->norm_bar_delta_x_prev = wnorm; 15770d8d27cSBarry Smith are set, but wnorm has not even been computed. 15870d8d27cSBarry Smith I don't know if this is the correct fix but by setting ynorm and wnorm to -1.0 at 15970d8d27cSBarry Smith least the linesearch object is kept in the state set by the SNESLineSearchReset_NLEQERR() call above 16070d8d27cSBarry Smith */ 16170d8d27cSBarry Smith ynorm = wnorm = -1.0; 162d4c6564cSPatrick Farrell break; 163d4c6564cSPatrick Farrell } 164d4c6564cSPatrick Farrell 165d4c6564cSPatrick Farrell /* Compute new trial iterate */ 1669566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -lambda, Y, X)); 1679566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); 168d4c6564cSPatrick Farrell 169d4c6564cSPatrick 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 */ 1709566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, G, W)); 1719566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 172d4c6564cSPatrick Farrell if (kspreason < 0) { 1739566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Solution for \\bar{delta x}^{k+1} failed.")); 174d4c6564cSPatrick Farrell } 175d4c6564cSPatrick Farrell 176d4c6564cSPatrick Farrell /* W now contains -bar_delta_x_curr. */ 177d4c6564cSPatrick Farrell 1789566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 179d4c6564cSPatrick Farrell if (monitor) { 1809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: norm of simplified Newton update: %14.12e\n", (double) wnorm)); 1829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 183d4c6564cSPatrick Farrell } 184d4c6564cSPatrick Farrell 185d4c6564cSPatrick Farrell /* compute the monitoring quantities theta and mudash. */ 186d4c6564cSPatrick Farrell 187d4c6564cSPatrick Farrell theta = wnorm / ynorm; 188d4c6564cSPatrick Farrell 1899566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -(1.0 - lambda), Y, W)); 1909566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); 191d4c6564cSPatrick Farrell 192d4c6564cSPatrick Farrell mudash = (0.5 * ynorm * lambda * lambda) / gnorm; 193d4c6564cSPatrick Farrell 194d4c6564cSPatrick Farrell /* Check for termination of the linesearch */ 195d4c6564cSPatrick Farrell if (theta >= 1.0) { 196d4c6564cSPatrick Farrell /* need to go around again with smaller lambda */ 197d4c6564cSPatrick Farrell if (monitor) { 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel)); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor," Line search: monotonicity check failed, ratio: %14.12e\n", (double) theta)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel)); 201d4c6564cSPatrick Farrell } 202d4c6564cSPatrick Farrell lambda = PetscMin(mudash, 0.5 * lambda); 203d4c6564cSPatrick Farrell lambda = PetscMax(lambda, minlambda); 204d4c6564cSPatrick Farrell /* continue through the loop, i.e. go back to regularity test */ 205bf8ddb9bSPatrick Farrell } else { 206d4c6564cSPatrick Farrell /* linesearch terminated */ 207d4c6564cSPatrick Farrell lambdadash = PetscMin(1.0, mudash); 208d4c6564cSPatrick Farrell 209d4c6564cSPatrick Farrell if (lambdadash == 1.0 && lambda == 1.0 && wnorm <= stol) { 210d4c6564cSPatrick Farrell /* store the updated state, X - Y - W, in G: 211d4c6564cSPatrick Farrell I need to keep W for the next linesearch */ 212*ef46b1a6SStefano Zampini PetscCall(VecWAXPY(G, -1.0, Y, X)); 2139566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, -1.0, W)); 214d4c6564cSPatrick Farrell break; 215d4c6564cSPatrick Farrell } 216d4c6564cSPatrick Farrell 217d4c6564cSPatrick Farrell /* Deuflhard suggests to add the following: 218d4c6564cSPatrick Farrell else if (lambdadash >= 4.0 * lambda) { 219d4c6564cSPatrick Farrell lambda = lambdadash; 220d4c6564cSPatrick Farrell } 221d4c6564cSPatrick Farrell to continue through the loop, i.e. go back to regularity test. 222d4c6564cSPatrick Farrell I deliberately exclude this, as I have practical experience of this 223d4c6564cSPatrick Farrell getting stuck in infinite loops (on e.g. an Allen--Cahn problem). */ 224d4c6564cSPatrick Farrell 225d4c6564cSPatrick Farrell else { 226d4c6564cSPatrick Farrell /* accept iterate without adding on, i.e. don't use bar_delta_x; 227d4c6564cSPatrick Farrell again, I need to keep W for the next linesearch */ 2289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -lambda, Y, X)); 229d4c6564cSPatrick Farrell break; 230d4c6564cSPatrick Farrell } 231d4c6564cSPatrick Farrell } 232d4c6564cSPatrick Farrell } 233d4c6564cSPatrick Farrell 234d4c6564cSPatrick Farrell if (linesearch->ops->viproject) { 2359566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->viproject)(snes, G)); 236d4c6564cSPatrick Farrell } 237d4c6564cSPatrick Farrell 238d4c6564cSPatrick Farrell /* W currently contains -bar_delta_u. Scale it so that it contains bar_delta_u. */ 2399566063dSJacob Faibussowitsch PetscCall(VecScale(W, -1.0)); 240d4c6564cSPatrick Farrell 241d4c6564cSPatrick Farrell /* postcheck */ 2429566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,G,&changed_y,&changed_w)); 243d4c6564cSPatrick Farrell if (changed_y || changed_w) { 2449566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_USER)); 2459566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Changing the search direction here doesn't make sense.\n")); 246d4c6564cSPatrick Farrell PetscFunctionReturn(0); 247d4c6564cSPatrick Farrell } 248d4c6564cSPatrick Farrell 249d4c6564cSPatrick Farrell /* copy the solution and information from this iteration over */ 250d4c6564cSPatrick Farrell nleqerr->norm_delta_x_prev = ynorm; 251d4c6564cSPatrick Farrell nleqerr->norm_bar_delta_x_prev = wnorm; 252d4c6564cSPatrick Farrell nleqerr->lambda_prev = lambda; 253d4c6564cSPatrick Farrell 2549566063dSJacob Faibussowitsch PetscCall(VecCopy(G, X)); 2559566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2569566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); 2579566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 2589566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetLambda(linesearch, lambda)); 2599566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, (ynorm < 0 ? PETSC_INFINITY : ynorm))); 260d4c6564cSPatrick Farrell PetscFunctionReturn(0); 261d4c6564cSPatrick Farrell } 262d4c6564cSPatrick Farrell 263d4c6564cSPatrick Farrell PetscErrorCode SNESLineSearchView_NLEQERR(SNESLineSearch linesearch, PetscViewer viewer) 264d4c6564cSPatrick Farrell { 265d4c6564cSPatrick Farrell PetscBool iascii; 266d4c6564cSPatrick Farrell SNESLineSearch_NLEQERR *nleqerr; 267d4c6564cSPatrick Farrell 268d4c6564cSPatrick Farrell PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 270d4c6564cSPatrick Farrell nleqerr = (SNESLineSearch_NLEQERR*)linesearch->data; 271d4c6564cSPatrick Farrell if (iascii) { 2729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " NLEQ-ERR affine-covariant linesearch")); 2739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " current local Lipschitz estimate omega=%e\n", (double)nleqerr->mu_curr)); 274d4c6564cSPatrick Farrell } 275d4c6564cSPatrick Farrell PetscFunctionReturn(0); 276d4c6564cSPatrick Farrell } 277d4c6564cSPatrick Farrell 278d4c6564cSPatrick Farrell static PetscErrorCode SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch) 279d4c6564cSPatrick Farrell { 280d4c6564cSPatrick Farrell PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(PetscFree(linesearch->data)); 282d4c6564cSPatrick Farrell PetscFunctionReturn(0); 283d4c6564cSPatrick Farrell } 284d4c6564cSPatrick Farrell 285d4c6564cSPatrick Farrell /*MC 286d4c6564cSPatrick Farrell SNESLINESEARCHNLEQERR - Error-oriented affine-covariant globalised Newton algorithm of Deuflhard (2011). 287d4c6564cSPatrick Farrell 288d4c6564cSPatrick Farrell This linesearch is intended for Newton-type methods which are affine covariant. Affine covariance 289d4c6564cSPatrick Farrell means that Newton's method will give the same iterations for F(x) = 0 and AF(x) = 0 for a nonsingular 290d4c6564cSPatrick Farrell matrix A. This is a fundamental property; the philosophy of this linesearch is that globalisations 291d4c6564cSPatrick Farrell of Newton's method should carefully preserve it. 292d4c6564cSPatrick Farrell 293d4c6564cSPatrick Farrell For a discussion of the theory behind this algorithm, see 294d4c6564cSPatrick Farrell 295d4c6564cSPatrick Farrell @book{deuflhard2011, 296d4c6564cSPatrick Farrell title={Newton Methods for Nonlinear Problems}, 297d4c6564cSPatrick Farrell author={Deuflhard, P.}, 298d4c6564cSPatrick Farrell volume={35}, 299d4c6564cSPatrick Farrell year={2011}, 300d4c6564cSPatrick Farrell publisher={Springer-Verlag}, 301d4c6564cSPatrick Farrell address={Berlin, Heidelberg} 302d4c6564cSPatrick Farrell } 303d4c6564cSPatrick Farrell 304d4c6564cSPatrick Farrell Pseudocode is given on page 148. 305d4c6564cSPatrick Farrell 306d4c6564cSPatrick Farrell Options Database Keys: 307d4c6564cSPatrick Farrell + -snes_linesearch_damping<1.0> - initial step length 308d4c6564cSPatrick Farrell - -snes_linesearch_minlambda<1e-12> - minimum step length allowed 309d4c6564cSPatrick Farrell 310d4c6564cSPatrick Farrell Contributed by Patrick Farrell <patrick.farrell@maths.ox.ac.uk> 311d4c6564cSPatrick Farrell 312d4c6564cSPatrick Farrell Level: advanced 313d4c6564cSPatrick Farrell 314d4c6564cSPatrick Farrell .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 315d4c6564cSPatrick Farrell M*/ 316d4c6564cSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch) 317d4c6564cSPatrick Farrell { 318d4c6564cSPatrick Farrell SNESLineSearch_NLEQERR *nleqerr; 319d4c6564cSPatrick Farrell 320d4c6564cSPatrick Farrell PetscFunctionBegin; 321d4c6564cSPatrick Farrell linesearch->ops->apply = SNESLineSearchApply_NLEQERR; 322d4c6564cSPatrick Farrell linesearch->ops->destroy = SNESLineSearchDestroy_NLEQERR; 323d4c6564cSPatrick Farrell linesearch->ops->setfromoptions = NULL; 324d4c6564cSPatrick Farrell linesearch->ops->reset = SNESLineSearchReset_NLEQERR; 325d4c6564cSPatrick Farrell linesearch->ops->view = SNESLineSearchView_NLEQERR; 326d4c6564cSPatrick Farrell linesearch->ops->setup = NULL; 327d4c6564cSPatrick Farrell 3289566063dSJacob Faibussowitsch PetscCall(PetscNewLog(linesearch,&nleqerr)); 329d4c6564cSPatrick Farrell 330d4c6564cSPatrick Farrell linesearch->data = (void*)nleqerr; 331d4c6564cSPatrick Farrell linesearch->max_its = 40; 3327cbffc34SPatrick Farrell SNESLineSearchReset_NLEQERR(linesearch); 333d4c6564cSPatrick Farrell PetscFunctionReturn(0); 334d4c6564cSPatrick Farrell } 335