Actual source code: linesearchbasic.c

petsc-3.3-p7 2013-05-11
  1: #include <petsc-private/linesearchimpl.h>
  2: #include <petsc-private/snesimpl.h>

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


 17:   SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);
 18:   SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);
 19:   SNESLineSearchGetLambda(linesearch, &lambda);
 20:   SNESLineSearchGetSNES(linesearch, &snes);
 21:   SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);

 23:   /* precheck */
 24:   SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);

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

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

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

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

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

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

 78:    Options Database Keys:
 79: +   -snes_linesearch_damping (1.0) damping parameter.
 80: -   -snes_linesearch_norms (true) whether to compute norms or not.

 82:    Notes:
 83:    For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
 84:    damping parameter may yield satisfactory but slow convergence despite
 85:    the simplicity of the line search.

 87:    Level: advanced

 89: .keywords: SNES, SNESLineSearch, damping

 91: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
 92: M*/
 94: {
 96:   linesearch->ops->apply          = SNESLineSearchApply_Basic;
 97:   linesearch->ops->destroy        = PETSC_NULL;
 98:   linesearch->ops->setfromoptions = PETSC_NULL;
 99:   linesearch->ops->reset          = PETSC_NULL;
100:   linesearch->ops->view           = PETSC_NULL;
101:   linesearch->ops->setup          = PETSC_NULL;
102:   return(0);
103: }