1bf7f4e0aSPeter Brune #include <private/linesearchimpl.h> 2bf7f4e0aSPeter Brune #include <private/snesimpl.h> 3bf7f4e0aSPeter Brune 4f40b411bSPeter Brune /*MC 56188f407SPeter Brune PetscLineSearchBasic - This routine is not a line search at all; 6bf7f4e0aSPeter Brune it simply uses the full step. Thus, this routine is intended 7bf7f4e0aSPeter Brune to serve as a template and is not recommended for general use. 8bf7f4e0aSPeter Brune 9bf7f4e0aSPeter Brune Level: advanced 10bf7f4e0aSPeter Brune 116188f407SPeter Brune .keywords: SNES, PetscLineSearch, damping 12bf7f4e0aSPeter Brune 136188f407SPeter Brune .seealso: PetscLineSearchCreate(), PetscLineSearchSetType() 14f40b411bSPeter Brune M*/ 15f40b411bSPeter Brune 16f40b411bSPeter Brune #undef __FUNCT__ 176188f407SPeter Brune #define __FUNCT__ "PetscLineSearchApply_Basic" 18f40b411bSPeter Brune 196188f407SPeter Brune PetscErrorCode PetscLineSearchApply_Basic(PetscLineSearch linesearch) 20bf7f4e0aSPeter Brune { 21bf7f4e0aSPeter Brune PetscBool changed_y, changed_w; 22bf7f4e0aSPeter Brune PetscErrorCode ierr; 236a388c36SPeter Brune Vec X, F, Y, W; 246a388c36SPeter Brune SNES snes; 256a388c36SPeter Brune PetscReal gnorm, xnorm, ynorm, lambda; 266a388c36SPeter Brune PetscBool domainerror; 27bf7f4e0aSPeter Brune 28bf7f4e0aSPeter Brune PetscFunctionBegin; 296a388c36SPeter Brune 306188f407SPeter Brune ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr); 316188f407SPeter Brune ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 326188f407SPeter Brune ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 336188f407SPeter Brune ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 346188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 35bf7f4e0aSPeter Brune 36bf7f4e0aSPeter Brune /* precheck */ 376188f407SPeter Brune ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 38bf7f4e0aSPeter Brune 39bf7f4e0aSPeter Brune /* update */ 406a388c36SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 41*9bd66eb0SPeter Brune if (linesearch->ops->viproject) { 42*9bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 43*9bd66eb0SPeter Brune } 44bf7f4e0aSPeter Brune 45bf7f4e0aSPeter Brune /* postcheck */ 466188f407SPeter Brune ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 47bf7f4e0aSPeter Brune if (changed_y) { 486a388c36SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 49*9bd66eb0SPeter Brune if (linesearch->ops->viproject) { 50*9bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 51*9bd66eb0SPeter Brune } 52bf7f4e0aSPeter Brune } 53bf7f4e0aSPeter Brune ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr); 546a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 556a388c36SPeter Brune if (domainerror) { 566188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE); 57bf7f4e0aSPeter Brune PetscFunctionReturn(0); 58bf7f4e0aSPeter Brune } 596a388c36SPeter Brune 60*9bd66eb0SPeter Brune ierr = VecNorm(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 61*9bd66eb0SPeter Brune ierr = VecNorm(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 62*9bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 63*9bd66eb0SPeter Brune linesearch->fnorm = gnorm; 64*9bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr); 65*9bd66eb0SPeter Brune } else { 66*9bd66eb0SPeter Brune ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr); 67*9bd66eb0SPeter Brune } 68*9bd66eb0SPeter Brune 69*9bd66eb0SPeter Brune /* 706188f407SPeter Brune ierr = PetscLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 71*9bd66eb0SPeter Brune */ 726a388c36SPeter Brune 73bf7f4e0aSPeter Brune /* copy the solution over */ 74bf7f4e0aSPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 75bf7f4e0aSPeter Brune PetscFunctionReturn(0); 76bf7f4e0aSPeter Brune } 77bf7f4e0aSPeter Brune 78bf7f4e0aSPeter Brune EXTERN_C_BEGIN 79bf7f4e0aSPeter Brune #undef __FUNCT__ 806188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate_Basic" 816188f407SPeter Brune PetscErrorCode PetscLineSearchCreate_Basic(PetscLineSearch linesearch) 82bf7f4e0aSPeter Brune { 83bf7f4e0aSPeter Brune PetscFunctionBegin; 846188f407SPeter Brune linesearch->ops->apply = PetscLineSearchApply_Basic; 85bf7f4e0aSPeter Brune linesearch->ops->destroy = PETSC_NULL; 86bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = PETSC_NULL; 87bf7f4e0aSPeter Brune linesearch->ops->reset = PETSC_NULL; 88bf7f4e0aSPeter Brune linesearch->ops->view = PETSC_NULL; 89bf7f4e0aSPeter Brune linesearch->ops->setup = PETSC_NULL; 90bf7f4e0aSPeter Brune PetscFunctionReturn(0); 91bf7f4e0aSPeter Brune } 92bf7f4e0aSPeter Brune EXTERN_C_END 93