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; 32bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr); 33b8ac21a2SPeter Brune fty_init = fty_old; 34b8ac21a2SPeter Brune 356a388c36SPeter Brune for (i = 0; i < max_its; i++) { 36bbd5d0b3SPeter Brune 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 } 43bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 44bbd5d0b3SPeter Brune 45bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 46bbd5d0b3SPeter Brune 47bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 48b8ac21a2SPeter Brune 49b8ac21a2SPeter Brune /* check for convergence */ 50b8ac21a2SPeter Brune if (PetscAbsReal(delLambda) < steptol*lambda) break; 51b8ac21a2SPeter Brune if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 52a7370c10SPeter Brune if (PetscAbsScalar(fty) < atol && i > 0) break; 536a388c36SPeter Brune if (monitor) { 546a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 55c69d1a72SBarry 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); 566a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 57bbd5d0b3SPeter Brune } 58bbd5d0b3SPeter Brune 59bbd5d0b3SPeter Brune /* compute the search direction */ 60b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { 61b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda; 62b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 63b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 64b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 65b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 66b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 67b8ac21a2SPeter Brune } 68b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 69b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 70b8ac21a2SPeter Brune s = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda; 71b8ac21a2SPeter Brune } else { 72b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 73b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 74b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 75b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 76b8ac21a2SPeter Brune } 77b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 78b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 79b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 80b8ac21a2SPeter Brune ierr = VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);CHKERRQ(ierr); 81b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 82b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 83b8ac21a2SPeter Brune } 84b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 85b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid2);CHKERRQ(ierr); 86b8ac21a2SPeter Brune s = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda); 87b8ac21a2SPeter Brune } 88b8ac21a2SPeter Brune /* if the solve is going in the wrong direction, fix it */ 89b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s; 90b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s); 91b8ac21a2SPeter Brune 92b8ac21a2SPeter Brune /* switch directions if we stepped out of bounds */ 93f5af7f23SKarl Rupp if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s); 94b8ac21a2SPeter Brune 95bbd5d0b3SPeter Brune if (PetscIsInfOrNanScalar(lambda_update)) break; 96f5af7f23SKarl Rupp if (lambda_update > maxstep) break; 97bbd5d0b3SPeter Brune 98bbd5d0b3SPeter Brune /* compute the new state of the line search */ 99bbd5d0b3SPeter Brune lambda_old = lambda; 100bbd5d0b3SPeter Brune lambda = lambda_update; 101bbd5d0b3SPeter Brune fty_old = fty; 102bbd5d0b3SPeter Brune } 103bbd5d0b3SPeter Brune /* construct the solution */ 104bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 105bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 1069bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1079bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1089bd66eb0SPeter Brune } 109bbd5d0b3SPeter Brune /* postcheck */ 1107b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 111bbd5d0b3SPeter Brune if (changed_y) { 112bbd5d0b3SPeter Brune ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 1139bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1149bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); 1159bd66eb0SPeter Brune } 116bbd5d0b3SPeter Brune } else { 117bbd5d0b3SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 118bbd5d0b3SPeter Brune } 119bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1206a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1216a388c36SPeter Brune if (domainerror) { 122f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 123bbd5d0b3SPeter Brune PetscFunctionReturn(0); 124bbd5d0b3SPeter Brune } 125bbd5d0b3SPeter Brune 126f1c6b773SPeter Brune ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 127f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 128bbd5d0b3SPeter Brune 129f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 1306a388c36SPeter Brune 1316a388c36SPeter Brune if (monitor) { 1326a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 133c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr); 1346a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 135bbd5d0b3SPeter Brune } 1366a388c36SPeter Brune if (lambda <= steptol) { 137f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 138bbd5d0b3SPeter Brune } 139bbd5d0b3SPeter Brune PetscFunctionReturn(0); 140bbd5d0b3SPeter Brune } 141bbd5d0b3SPeter Brune 142bbd5d0b3SPeter Brune #undef __FUNCT__ 143f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_CP" 144954494b2SJed Brown /*MC 1451a4f838cSPeter Brune SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some 146cd7522eaSPeter Brune artificial G(x) for which the SNESFunction F(x) = grad G(x). Therefore, this line search seeks 147db609ea7SPeter Brune to find roots of dot(F, Y) via a secant method. 148cd7522eaSPeter Brune 149cd7522eaSPeter Brune Options Database Keys: 150db609ea7SPeter Brune + -snes_linesearch_minlambda - the minimum acceptable lambda 151cd7522eaSPeter Brune . -snes_linesearch_damping - initial trial step length 152db609ea7SPeter Brune - -snes_linesearch_max_it - the maximum number of secant steps performed. 153cd7522eaSPeter Brune 154cd7522eaSPeter Brune Notes: 155cd7522eaSPeter Brune This method is the preferred line search for SNESQN and SNESNCG. 156954494b2SJed Brown 157954494b2SJed Brown Level: advanced 158954494b2SJed Brown 159f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 160954494b2SJed Brown 161f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 162954494b2SJed Brown M*/ 163*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 164bbd5d0b3SPeter Brune { 165bbd5d0b3SPeter Brune PetscFunctionBegin; 166f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP; 1670298fd71SBarry Smith linesearch->ops->destroy = NULL; 1680298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1690298fd71SBarry Smith linesearch->ops->reset = NULL; 1700298fd71SBarry Smith linesearch->ops->view = NULL; 1710298fd71SBarry Smith linesearch->ops->setup = NULL; 172b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_LINEAR; 173a491bab8SPeter Brune 174a491bab8SPeter Brune linesearch->max_its = 1; 175bbd5d0b3SPeter Brune PetscFunctionReturn(0); 176bbd5d0b3SPeter Brune } 177