xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision 954494b2f105144292715050ea1b8966feeb7ce6)
1bf7f4e0aSPeter Brune #include <private/linesearchimpl.h>
2bf7f4e0aSPeter Brune #include <private/snesimpl.h>
3bf7f4e0aSPeter Brune 
4f40b411bSPeter Brune #undef __FUNCT__
56188f407SPeter Brune #define __FUNCT__ "PetscLineSearchApply_Basic"
6*954494b2SJed Brown static PetscErrorCode  PetscLineSearchApply_Basic(PetscLineSearch linesearch)
7bf7f4e0aSPeter Brune {
8bf7f4e0aSPeter Brune   PetscBool      changed_y, changed_w;
9bf7f4e0aSPeter Brune   PetscErrorCode ierr;
106a388c36SPeter Brune   Vec            X, F, Y, W;
116a388c36SPeter Brune   SNES           snes;
126a388c36SPeter Brune   PetscReal      gnorm, xnorm, ynorm, lambda;
136a388c36SPeter Brune   PetscBool      domainerror;
14bf7f4e0aSPeter Brune 
15bf7f4e0aSPeter Brune   PetscFunctionBegin;
166a388c36SPeter Brune 
176188f407SPeter Brune   ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr);
186188f407SPeter Brune   ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
196188f407SPeter Brune   ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
206188f407SPeter Brune   ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
216188f407SPeter Brune   ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
22bf7f4e0aSPeter Brune 
23bf7f4e0aSPeter Brune   /* precheck */
246188f407SPeter Brune   ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
25bf7f4e0aSPeter Brune 
26bf7f4e0aSPeter Brune   /* update */
276a388c36SPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
289bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
299bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
309bd66eb0SPeter Brune   }
31bf7f4e0aSPeter Brune 
32bf7f4e0aSPeter Brune   /* postcheck */
336188f407SPeter Brune   ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
34bf7f4e0aSPeter Brune   if (changed_y) {
356a388c36SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
369bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
379bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
389bd66eb0SPeter Brune     }
39bf7f4e0aSPeter Brune   }
40bf7f4e0aSPeter Brune   ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
416a388c36SPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
426a388c36SPeter Brune   if (domainerror) {
436188f407SPeter Brune     ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);
44bf7f4e0aSPeter Brune     PetscFunctionReturn(0);
45bf7f4e0aSPeter Brune   }
466a388c36SPeter Brune 
479bd66eb0SPeter Brune   ierr = VecNorm(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
489bd66eb0SPeter Brune   ierr = VecNorm(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
499bd66eb0SPeter Brune   if (linesearch->ops->vinorm) {
509bd66eb0SPeter Brune     linesearch->fnorm = gnorm;
519bd66eb0SPeter Brune     ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr);
529bd66eb0SPeter Brune   } else {
539bd66eb0SPeter Brune     ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr);
549bd66eb0SPeter Brune   }
559bd66eb0SPeter Brune 
569bd66eb0SPeter Brune   /*
576188f407SPeter Brune   ierr = PetscLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
589bd66eb0SPeter Brune    */
596a388c36SPeter Brune 
60bf7f4e0aSPeter Brune   /* copy the solution over */
61bf7f4e0aSPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
62bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
63bf7f4e0aSPeter Brune }
64bf7f4e0aSPeter Brune 
65bf7f4e0aSPeter Brune #undef __FUNCT__
666188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate_Basic"
67*954494b2SJed Brown /*MC
68*954494b2SJed Brown    PETSCLINESEARCHBASIC - This routine is not a line search at all;
69*954494b2SJed Brown    it simply uses the full step.  Thus, this routine is intended
70*954494b2SJed Brown    to serve as a template and is not recommended for general use.
71*954494b2SJed Brown 
72*954494b2SJed Brown    Level: advanced
73*954494b2SJed Brown 
74*954494b2SJed Brown .keywords: SNES, PetscLineSearch, damping
75*954494b2SJed Brown 
76*954494b2SJed Brown .seealso: PetscLineSearchCreate(), PetscLineSearchSetType()
77*954494b2SJed Brown M*/
78*954494b2SJed Brown PETSC_EXTERN_C PetscErrorCode PetscLineSearchCreate_Basic(PetscLineSearch linesearch)
79bf7f4e0aSPeter Brune {
80bf7f4e0aSPeter Brune   PetscFunctionBegin;
816188f407SPeter Brune   linesearch->ops->apply          = PetscLineSearchApply_Basic;
82bf7f4e0aSPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
83bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
84bf7f4e0aSPeter Brune   linesearch->ops->reset          = PETSC_NULL;
85bf7f4e0aSPeter Brune   linesearch->ops->view           = PETSC_NULL;
86bf7f4e0aSPeter Brune   linesearch->ops->setup          = PETSC_NULL;
87bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
88bf7f4e0aSPeter Brune }
89