1*b45d2f2cSJed 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; 216a388c36SPeter Brune 22f1c6b773SPeter Brune ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr); 23f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 24f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 25f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 26f1c6b773SPeter Brune ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its);CHKERRQ(ierr); 27f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 28f1c6b773SPeter Brune ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 296a388c36SPeter Brune 30bbd5d0b3SPeter Brune /* precheck */ 31f1c6b773SPeter Brune ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 32bbd5d0b3SPeter Brune lambda_old = 0.0; 33bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr); 34b8ac21a2SPeter Brune fty_init = fty_old; 35b8ac21a2SPeter Brune 366a388c36SPeter Brune for (i = 0; i < max_its; i++) { 37bbd5d0b3SPeter Brune 38bbd5d0b3SPeter Brune /* compute the norm at lambda */ 39bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 40bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 419bd66eb0SPeter Brune if (linesearch->ops->viproject) { 429bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 439bd66eb0SPeter Brune } 44bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 45bbd5d0b3SPeter Brune 46bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 47bbd5d0b3SPeter Brune 48bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 49b8ac21a2SPeter Brune 50b8ac21a2SPeter Brune /* check for convergence */ 51b8ac21a2SPeter Brune if (PetscAbsReal(delLambda) < steptol*lambda) break; 52b8ac21a2SPeter Brune if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 53a7370c10SPeter Brune if (PetscAbsScalar(fty) < atol && i > 0) break; 546a388c36SPeter Brune if (monitor) { 556a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 566a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g], ftys = [%g, %g]\n", 57bbd5d0b3SPeter Brune lambda, lambda_old, PetscRealPart(fty), PetscRealPart(fty_old));CHKERRQ(ierr); 586a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 59bbd5d0b3SPeter Brune } 60bbd5d0b3SPeter Brune 61bbd5d0b3SPeter Brune /* compute the search direction */ 62b8ac21a2SPeter Brune if (linesearch->order == SNES_LINESEARCH_LINEAR) { 63b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda; 64b8ac21a2SPeter Brune } else if (linesearch->order == SNES_LINESEARCH_QUADRATIC) { 65b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 66b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 67b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 68b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 69b8ac21a2SPeter Brune } 70b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 71b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 72b8ac21a2SPeter Brune s = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda; 73b8ac21a2SPeter Brune } else { 74b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 75b8ac21a2SPeter Brune ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 76b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 77b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 78b8ac21a2SPeter Brune } 79b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 80b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 81b8ac21a2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 82b8ac21a2SPeter Brune ierr = VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);CHKERRQ(ierr); 83b8ac21a2SPeter Brune if (linesearch->ops->viproject) { 84b8ac21a2SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 85b8ac21a2SPeter Brune } 86b8ac21a2SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 87b8ac21a2SPeter Brune ierr = VecDot(F, Y, &fty_mid2);CHKERRQ(ierr); 88b8ac21a2SPeter Brune s = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda); 89b8ac21a2SPeter Brune } 90b8ac21a2SPeter Brune /* if the solve is going in the wrong direction, fix it */ 91b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s; 92b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s); 93b8ac21a2SPeter Brune 94b8ac21a2SPeter Brune /* switch directions if we stepped out of bounds */ 95b8ac21a2SPeter Brune if (lambda_update < steptol) { 96b8ac21a2SPeter Brune lambda_update = lambda + PetscRealPart(fty / s); 97b8ac21a2SPeter Brune } 98b8ac21a2SPeter Brune 99bbd5d0b3SPeter Brune if (PetscIsInfOrNanScalar(lambda_update)) break; 1006a388c36SPeter Brune if (lambda_update > maxstep) { 101bbd5d0b3SPeter Brune break; 102bbd5d0b3SPeter Brune } 103bbd5d0b3SPeter Brune 104bbd5d0b3SPeter Brune /* compute the new state of the line search */ 105bbd5d0b3SPeter Brune lambda_old = lambda; 106bbd5d0b3SPeter Brune lambda = lambda_update; 107bbd5d0b3SPeter Brune fty_old = fty; 108bbd5d0b3SPeter Brune } 109bbd5d0b3SPeter Brune /* construct the solution */ 110bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 111bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 1129bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1139bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1149bd66eb0SPeter Brune } 115bbd5d0b3SPeter Brune /* postcheck */ 116f1c6b773SPeter Brune ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 117bbd5d0b3SPeter Brune if (changed_y) { 118bbd5d0b3SPeter Brune ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 1199bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1209bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); 1219bd66eb0SPeter Brune } 122bbd5d0b3SPeter Brune } else { 123bbd5d0b3SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 124bbd5d0b3SPeter Brune } 125bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1266a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1276a388c36SPeter Brune if (domainerror) { 128f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 129bbd5d0b3SPeter Brune PetscFunctionReturn(0); 130bbd5d0b3SPeter Brune } 131bbd5d0b3SPeter Brune 132f1c6b773SPeter Brune ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 133f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 134bbd5d0b3SPeter Brune 135f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 1366a388c36SPeter Brune 1376a388c36SPeter Brune if (monitor) { 1386a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 1396a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", lambda, gnorm);CHKERRQ(ierr); 1406a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 141bbd5d0b3SPeter Brune } 1426a388c36SPeter Brune if (lambda <= steptol) { 143f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 144bbd5d0b3SPeter Brune } 145bbd5d0b3SPeter Brune PetscFunctionReturn(0); 146bbd5d0b3SPeter Brune } 147bbd5d0b3SPeter Brune 148bbd5d0b3SPeter Brune #undef __FUNCT__ 149f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_CP" 150954494b2SJed Brown /*MC 151cd7522eaSPeter Brune SNES_LINESEARCH_CP - Critical point line search. This line search assumes that there exists some 152cd7522eaSPeter Brune artificial G(x) for which the SNESFunction F(x) = grad G(x). Therefore, this line search seeks 153cd7522eaSPeter Brune to find roots of f^ty via a secant method. 154cd7522eaSPeter Brune 155cd7522eaSPeter Brune Options Database Keys: 156cd7522eaSPeter Brune . -snes_linesearch_damping - initial trial step length 157cd7522eaSPeter Brune 158cd7522eaSPeter Brune Notes: 159cd7522eaSPeter Brune This method is the preferred line search for SNESQN and SNESNCG. 160954494b2SJed Brown 161954494b2SJed Brown Level: advanced 162954494b2SJed Brown 163f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 164954494b2SJed Brown 165f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 166954494b2SJed Brown M*/ 167f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 168bbd5d0b3SPeter Brune { 169bbd5d0b3SPeter Brune PetscFunctionBegin; 170f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP; 171bbd5d0b3SPeter Brune linesearch->ops->destroy = PETSC_NULL; 172bbd5d0b3SPeter Brune linesearch->ops->setfromoptions = PETSC_NULL; 173bbd5d0b3SPeter Brune linesearch->ops->reset = PETSC_NULL; 174bbd5d0b3SPeter Brune linesearch->ops->view = PETSC_NULL; 175bbd5d0b3SPeter Brune linesearch->ops->setup = PETSC_NULL; 176b8ac21a2SPeter Brune linesearch->order = SNES_LINESEARCH_LINEAR; 177bbd5d0b3SPeter Brune PetscFunctionReturn(0); 178bbd5d0b3SPeter Brune } 179