1 #include <petsc-private/linesearchimpl.h> 2 #include <petsc-private/snesimpl.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESLineSearchApply_Basic" 6 static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch) 7 { 8 PetscBool changed_y, changed_w; 9 PetscErrorCode ierr; 10 Vec X, F, Y, W; 11 SNES snes; 12 PetscReal gnorm, xnorm, ynorm, lambda; 13 PetscBool domainerror; 14 15 PetscFunctionBegin; 16 ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);CHKERRQ(ierr); 17 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 18 ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 19 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 20 ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 21 22 /* precheck */ 23 ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 24 25 /* update */ 26 ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 27 if (linesearch->ops->viproject) { 28 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 29 } 30 31 /* postcheck */ 32 ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 33 if (changed_y) { 34 ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 35 if (linesearch->ops->viproject) { 36 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 37 } 38 } 39 if (linesearch->norms || snes->iter < snes->max_its-1) { 40 ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr); 41 ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 42 if (domainerror) { 43 ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 44 PetscFunctionReturn(0); 45 } 46 } 47 48 if (linesearch->norms) { 49 if (!linesearch->ops->vinorm) VecNormBegin(F, NORM_2, &linesearch->fnorm); 50 ierr = VecNormBegin(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 51 ierr = VecNormBegin(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 52 if (!linesearch->ops->vinorm) VecNormEnd(F, NORM_2, &linesearch->fnorm); 53 ierr = VecNormEnd(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 54 ierr = VecNormEnd(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 55 56 if (linesearch->ops->vinorm) { 57 linesearch->fnorm = gnorm; 58 59 ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr); 60 } else { 61 ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr); 62 } 63 } 64 65 /* copy the solution over */ 66 ierr = VecCopy(W, X);CHKERRQ(ierr); 67 PetscFunctionReturn(0); 68 } 69 70 #undef __FUNCT__ 71 #define __FUNCT__ "SNESLineSearchCreate_Basic" 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 (SNESNEWTONLS), on 76 well-behaved problems. 77 78 Options Database Keys: 79 + -snes_linesearch_damping <damping> search vector is scaled by this amount, default is 1.0 80 - -snes_linesearch_norms <flag> whether to compute norms or not, default is true 81 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. 86 87 Level: advanced 88 89 .keywords: SNES, SNESLineSearch, damping 90 91 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 92 M*/ 93 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) 94 { 95 PetscFunctionBegin; 96 linesearch->ops->apply = SNESLineSearchApply_Basic; 97 linesearch->ops->destroy = NULL; 98 linesearch->ops->setfromoptions = NULL; 99 linesearch->ops->reset = NULL; 100 linesearch->ops->view = NULL; 101 linesearch->ops->setup = NULL; 102 PetscFunctionReturn(0); 103 } 104