1b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h> 26a388c36SPeter Brune #include <petscsnes.h> 3bbd5d0b3SPeter Brune 4bbd5d0b3SPeter Brune #undef __FUNCT__ 5f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_CP" 6f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch) 7bbd5d0b3SPeter Brune { 86a388c36SPeter Brune PetscBool changed_y, changed_w, domainerror; 9bbd5d0b3SPeter Brune PetscErrorCode ierr; 106a388c36SPeter Brune Vec X, Y, F, W; 116a388c36SPeter Brune SNES snes; 12516fe3c3SPeter Brune PetscReal xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep; 13bbd5d0b3SPeter Brune 14bbd5d0b3SPeter Brune PetscReal lambda, lambda_old, lambda_update, delLambda; 15b8ac21a2SPeter Brune PetscScalar fty, fty_init, fty_old, fty_mid1, fty_mid2, s; 166a388c36SPeter Brune PetscInt i, max_its; 176a388c36SPeter Brune 186a388c36SPeter Brune PetscViewer monitor; 19bbd5d0b3SPeter Brune 20bbd5d0b3SPeter Brune PetscFunctionBegin; 210298fd71SBarry Smith ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);CHKERRQ(ierr); 22f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 23f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 24f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 25f1c6b773SPeter Brune ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its);CHKERRQ(ierr); 26f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 27f1c6b773SPeter Brune ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 286a388c36SPeter Brune 29bbd5d0b3SPeter Brune /* precheck */ 307b1df9c1SPeter Brune ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 31bbd5d0b3SPeter Brune lambda_old = 0.0; 3262893cf2SPeter Brune 33e0e0c2afSPeter Brune ierr = VecDot(F,Y,&fty_old);CHKERRQ(ierr); 3462893cf2SPeter Brune fty_init = fty_old; 35*c0b2db79SJed Brown 36*c0b2db79SJed Brown for (i = 0; i < max_its; i++) { 37bbd5d0b3SPeter Brune /* compute the norm at lambda */ 38bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 39bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 409bd66eb0SPeter Brune if (linesearch->ops->viproject) { 419bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 429bd66eb0SPeter Brune } 43fb8e56e0SPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr); 44e0e0c2afSPeter Brune ierr = VecDot(F,Y,&fty);CHKERRQ(ierr); 45bbd5d0b3SPeter Brune 46bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 47b8ac21a2SPeter Brune 48b8ac21a2SPeter Brune /* check for convergence */ 49b8ac21a2SPeter Brune if (PetscAbsReal(delLambda) < steptol*lambda) break; 50b8ac21a2SPeter Brune if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 51a7370c10SPeter Brune if (PetscAbsScalar(fty) < atol && i > 0) break; 526a388c36SPeter Brune if (monitor) { 536a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 54c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old));CHKERRQ(ierr); 556a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 56bbd5d0b3SPeter Brune } 57bbd5d0b3SPeter Brune 58bbd5d0b3SPeter Brune /* compute the search direction */ 59b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { 60b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda; 61b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 62b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 63b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 64b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 65b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 66b8ac21a2SPeter Brune } 67fb8e56e0SPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr); 68b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 69b8ac21a2SPeter Brune s = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda; 70b8ac21a2SPeter Brune } else { 71b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 72b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 73b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 74b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 75b8ac21a2SPeter Brune } 76fb8e56e0SPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr); 77e0e0c2afSPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 78b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 79b8ac21a2SPeter Brune ierr = VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);CHKERRQ(ierr); 80b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 81b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 82b8ac21a2SPeter Brune } 83fb8e56e0SPeter Brune ierr = (*linesearch->ops->snesfunc)(snes, W, F);CHKERRQ(ierr); 84e0e0c2afSPeter Brune ierr = VecDot(F, Y, &fty_mid2);CHKERRQ(ierr); 85b8ac21a2SPeter Brune s = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda); 86b8ac21a2SPeter Brune } 87b8ac21a2SPeter Brune /* if the solve is going in the wrong direction, fix it */ 88b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s; 89b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s); 90b8ac21a2SPeter Brune 91b8ac21a2SPeter Brune /* switch directions if we stepped out of bounds */ 92f5af7f23SKarl Rupp if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s); 93b8ac21a2SPeter Brune 94bbd5d0b3SPeter Brune if (PetscIsInfOrNanScalar(lambda_update)) break; 95f5af7f23SKarl Rupp if (lambda_update > maxstep) break; 96bbd5d0b3SPeter Brune 97bbd5d0b3SPeter Brune /* compute the new state of the line search */ 98bbd5d0b3SPeter Brune lambda_old = lambda; 99bbd5d0b3SPeter Brune lambda = lambda_update; 100bbd5d0b3SPeter Brune fty_old = fty; 101bbd5d0b3SPeter Brune } 102bbd5d0b3SPeter Brune /* construct the solution */ 103bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 104bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 1059bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1069bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1079bd66eb0SPeter Brune } 108bbd5d0b3SPeter Brune /* postcheck */ 1097b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 110bbd5d0b3SPeter Brune if (changed_y) { 111bbd5d0b3SPeter Brune ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 1129bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1139bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); 1149bd66eb0SPeter Brune } 115bbd5d0b3SPeter Brune } else { 116bbd5d0b3SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 117bbd5d0b3SPeter Brune } 118fb8e56e0SPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,X,F);CHKERRQ(ierr); 1196a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1206a388c36SPeter Brune if (domainerror) { 121f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 122bbd5d0b3SPeter Brune PetscFunctionReturn(0); 123bbd5d0b3SPeter Brune } 124bbd5d0b3SPeter Brune 125f1c6b773SPeter Brune ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 126f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 127bbd5d0b3SPeter Brune 128f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 1296a388c36SPeter Brune 1306a388c36SPeter Brune if (monitor) { 1316a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 132c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr); 1336a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 134bbd5d0b3SPeter Brune } 1356a388c36SPeter Brune if (lambda <= steptol) { 136f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 137bbd5d0b3SPeter Brune } 138bbd5d0b3SPeter Brune PetscFunctionReturn(0); 139bbd5d0b3SPeter Brune } 140bbd5d0b3SPeter Brune 141bbd5d0b3SPeter Brune #undef __FUNCT__ 142f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_CP" 143954494b2SJed Brown /*MC 1441a4f838cSPeter Brune SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some 145cd7522eaSPeter Brune artificial G(x) for which the SNESFunction F(x) = grad G(x). Therefore, this line search seeks 146db609ea7SPeter Brune to find roots of dot(F, Y) via a secant method. 147cd7522eaSPeter Brune 148cd7522eaSPeter Brune Options Database Keys: 149db609ea7SPeter Brune + -snes_linesearch_minlambda - the minimum acceptable lambda 150cd7522eaSPeter Brune . -snes_linesearch_damping - initial trial step length 151db609ea7SPeter Brune - -snes_linesearch_max_it - the maximum number of secant steps performed. 152cd7522eaSPeter Brune 153cd7522eaSPeter Brune Notes: 154cd7522eaSPeter Brune This method is the preferred line search for SNESQN and SNESNCG. 155954494b2SJed Brown 156954494b2SJed Brown Level: advanced 157954494b2SJed Brown 158f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 159954494b2SJed Brown 160f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 161954494b2SJed Brown M*/ 1628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 163bbd5d0b3SPeter Brune { 164bbd5d0b3SPeter Brune PetscFunctionBegin; 165f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP; 1660298fd71SBarry Smith linesearch->ops->destroy = NULL; 1670298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1680298fd71SBarry Smith linesearch->ops->reset = NULL; 1690298fd71SBarry Smith linesearch->ops->view = NULL; 1700298fd71SBarry Smith linesearch->ops->setup = NULL; 171b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_LINEAR; 172a491bab8SPeter Brune 173a491bab8SPeter Brune linesearch->max_its = 1; 174bbd5d0b3SPeter Brune PetscFunctionReturn(0); 175bbd5d0b3SPeter Brune } 176