Actual source code: linesearchbasic.c

  1: #include <petsc/private/linesearchimpl.h>
  2: #include <petsc/private/snesimpl.h>

  4: static PetscErrorCode  SNESLineSearchApply_Basic(SNESLineSearch linesearch)
  5: {
  6:   PetscBool      changed_y, changed_w;
  7:   Vec            X, F, Y, W;
  8:   SNES           snes;
  9:   PetscReal      gnorm, xnorm, ynorm, lambda;
 10:   PetscBool      domainerror;

 12:   SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);
 13:   SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);
 14:   SNESLineSearchGetLambda(linesearch, &lambda);
 15:   SNESLineSearchGetSNES(linesearch, &snes);
 16:   SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);

 18:   /* precheck */
 19:   SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);

 21:   /* update */
 22:   VecWAXPY(W,-lambda,Y,X);
 23:   if (linesearch->ops->viproject) {
 24:     (*linesearch->ops->viproject)(snes, W);
 25:   }

 27:   /* postcheck */
 28:   SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);
 29:   if (changed_y) {
 30:     VecWAXPY(W,-lambda,Y,X);
 31:     if (linesearch->ops->viproject) {
 32:       (*linesearch->ops->viproject)(snes, W);
 33:     }
 34:   }
 35:   if (linesearch->norms || snes->iter < snes->max_its-1) {
 36:     (*linesearch->ops->snesfunc)(snes,W,F);
 37:     SNESGetFunctionDomainError(snes, &domainerror);
 38:     if (domainerror) {
 39:       SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN);
 40:       return 0;
 41:     }
 42:   }

 44:   if (linesearch->norms) {
 45:     if (!linesearch->ops->vinorm) VecNormBegin(F, NORM_2, &linesearch->fnorm);
 46:     VecNormBegin(Y, NORM_2, &linesearch->ynorm);
 47:     VecNormBegin(W, NORM_2, &linesearch->xnorm);
 48:     if (!linesearch->ops->vinorm) VecNormEnd(F, NORM_2, &linesearch->fnorm);
 49:     VecNormEnd(Y, NORM_2, &linesearch->ynorm);
 50:     VecNormEnd(W, NORM_2, &linesearch->xnorm);

 52:     if (linesearch->ops->vinorm) {
 53:       linesearch->fnorm = gnorm;

 55:       (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);
 56:     } else {
 57:       VecNorm(F,NORM_2,&linesearch->fnorm);
 58:     }
 59:   }

 61:   /* copy the solution over */
 62:   VecCopy(W, X);
 63:   return 0;
 64: }

 66: /*MC
 67:    SNESLINESEARCHBASIC - This line search implementation is not a line
 68:    search at all; it simply uses the full step.  Thus, this routine is intended
 69:    for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on
 70:    well-behaved problems.

 72:    Options Database Keys:
 73: +   -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0
 74: -   -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms())

 76:    Notes:
 77:    For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
 78:    damping parameter may yield satisfactory but slow convergence despite
 79:    the simplicity of the line search.

 81:    Level: advanced

 83: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType(), SNESLineSearchSetDamping(), SNESLineSearchSetComputeNorms()
 84: M*/
 85: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
 86: {
 87:   linesearch->ops->apply          = SNESLineSearchApply_Basic;
 88:   linesearch->ops->destroy        = NULL;
 89:   linesearch->ops->setfromoptions = NULL;
 90:   linesearch->ops->reset          = NULL;
 91:   linesearch->ops->view           = NULL;
 92:   linesearch->ops->setup          = NULL;
 93:   return 0;
 94: }