Actual source code: linesearchbasic.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  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;
  8:   Vec            X, F, Y, W;
  9:   SNES           snes;
 10:   PetscReal      gnorm, xnorm, ynorm, lambda;
 11:   PetscBool      domainerror;

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

 20:   /* precheck */
 21:   SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);

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

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

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

 54:     if (linesearch->ops->vinorm) {
 55:       linesearch->fnorm = gnorm;

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

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

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

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

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

 83:    Level: advanced

 85: .keywords: SNES, SNESLineSearch, damping

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