1bbd5d0b3SPeter Brune #include <private/linesearchimpl.h> 2*6a388c36SPeter Brune #include <petscsnes.h> 3bbd5d0b3SPeter Brune 4bbd5d0b3SPeter Brune #undef __FUNCT__ 5bbd5d0b3SPeter Brune #define __FUNCT__ "LineSearchApply_CP" 6bbd5d0b3SPeter Brune 7bbd5d0b3SPeter Brune /*@C 8bbd5d0b3SPeter Brune LineSearchCP - This routine is not a line search at all; 9bbd5d0b3SPeter Brune it simply uses the full step. Thus, this routine is intended 10bbd5d0b3SPeter Brune to serve as a template and is not recommended for general use. 11bbd5d0b3SPeter Brune 12bbd5d0b3SPeter Brune Logically Collective on SNES and Vec 13bbd5d0b3SPeter Brune 14bbd5d0b3SPeter Brune Input Parameters: 15bbd5d0b3SPeter Brune + snes - nonlinear context 16bbd5d0b3SPeter Brune . lsctx - optional context for line search (not used here) 17bbd5d0b3SPeter Brune . x - current iterate 18bbd5d0b3SPeter Brune . f - residual evaluated at x 19bbd5d0b3SPeter Brune . y - search direction 20bbd5d0b3SPeter Brune . fnorm - 2-norm of f 21bbd5d0b3SPeter Brune - xnorm - norm of x if known, otherwise 0 22bbd5d0b3SPeter Brune 23bbd5d0b3SPeter Brune Output Parameters: 24bbd5d0b3SPeter Brune + g - residual evaluated at new iterate y 25bbd5d0b3SPeter Brune . w - new iterate 26bbd5d0b3SPeter Brune . gnorm - 2-norm of g 27bbd5d0b3SPeter Brune . ynorm - 2-norm of search length 28bbd5d0b3SPeter Brune - flag - PETSC_TRUE on success, PETSC_FALSE on failure 29bbd5d0b3SPeter Brune 30bbd5d0b3SPeter Brune Options Database Key: 31bbd5d0b3SPeter Brune . -snes_ls basic - Activates SNESLineSearchNo() 32bbd5d0b3SPeter Brune 33bbd5d0b3SPeter Brune Level: advanced 34bbd5d0b3SPeter Brune 35bbd5d0b3SPeter Brune .keywords: SNES, nonlinear, line search, cubic 36bbd5d0b3SPeter Brune 37bbd5d0b3SPeter Brune .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 38bbd5d0b3SPeter Brune SNESLineSearchSet(), SNESLineSearchNoNorms() 39bbd5d0b3SPeter Brune @*/ 40bbd5d0b3SPeter Brune PetscErrorCode LineSearchApply_CP(LineSearch linesearch) 41bbd5d0b3SPeter Brune { 42bbd5d0b3SPeter Brune 43*6a388c36SPeter Brune PetscBool changed_y, changed_w, domainerror; 44bbd5d0b3SPeter Brune PetscErrorCode ierr; 45*6a388c36SPeter Brune Vec X, Y, F, W; 46*6a388c36SPeter Brune SNES snes; 47*6a388c36SPeter Brune PetscReal xnorm, ynorm, gnorm, steptol, maxstep; 48bbd5d0b3SPeter Brune 49bbd5d0b3SPeter Brune PetscReal lambda, lambda_old, lambda_update, delLambda; 50bbd5d0b3SPeter Brune PetscScalar fty, fty_old; 51*6a388c36SPeter Brune PetscInt i, max_its; 52*6a388c36SPeter Brune 53*6a388c36SPeter Brune PetscViewer monitor; 54bbd5d0b3SPeter Brune 55bbd5d0b3SPeter Brune PetscFunctionBegin; 56*6a388c36SPeter Brune 57*6a388c36SPeter Brune ierr = LineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr); 58*6a388c36SPeter Brune ierr = LineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 59*6a388c36SPeter Brune ierr = LineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 60*6a388c36SPeter Brune ierr = LineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 61*6a388c36SPeter Brune ierr = LineSearchGetMaxIts(linesearch, &max_its);CHKERRQ(ierr); 62*6a388c36SPeter Brune ierr = LineSearchGetStepTolerance(linesearch, &steptol);CHKERRQ(ierr); 63*6a388c36SPeter Brune ierr = LineSearchGetMaxStep(linesearch, &maxstep);CHKERRQ(ierr); 64*6a388c36SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 65*6a388c36SPeter Brune ierr = LineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 66*6a388c36SPeter Brune 67bbd5d0b3SPeter Brune /* precheck */ 68bbd5d0b3SPeter Brune ierr = LineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 69bbd5d0b3SPeter Brune lambda_old = 0.0; 70bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr); 71*6a388c36SPeter Brune for (i = 0; i < max_its; i++) { 72bbd5d0b3SPeter Brune 73bbd5d0b3SPeter Brune /* compute the norm at lambda */ 74bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 75bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 76bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 77bbd5d0b3SPeter Brune 78bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 79bbd5d0b3SPeter Brune 80bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 81*6a388c36SPeter Brune if (PetscAbsReal(delLambda) < steptol) break; 82*6a388c36SPeter Brune if (monitor) { 83*6a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 84*6a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g], ftys = [%g, %g]\n", 85bbd5d0b3SPeter Brune lambda, lambda_old, PetscRealPart(fty), PetscRealPart(fty_old));CHKERRQ(ierr); 86*6a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 87bbd5d0b3SPeter Brune } 88bbd5d0b3SPeter Brune 89bbd5d0b3SPeter Brune /* compute the search direction */ 90bbd5d0b3SPeter Brune lambda_update = PetscRealPart((fty*lambda_old - fty_old*lambda) / (fty - fty_old)); 91bbd5d0b3SPeter Brune if (PetscIsInfOrNanScalar(lambda_update)) break; 92*6a388c36SPeter Brune if (lambda_update > maxstep) { 93bbd5d0b3SPeter Brune break; 94bbd5d0b3SPeter Brune } 95bbd5d0b3SPeter Brune 96bbd5d0b3SPeter Brune /* compute the new state of the line search */ 97bbd5d0b3SPeter Brune lambda_old = lambda; 98bbd5d0b3SPeter Brune lambda = lambda_update; 99bbd5d0b3SPeter Brune fty_old = fty; 100bbd5d0b3SPeter Brune } 101bbd5d0b3SPeter Brune /* construct the solution */ 102bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 103bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 104bbd5d0b3SPeter Brune 105bbd5d0b3SPeter Brune /* postcheck */ 106bbd5d0b3SPeter Brune ierr = LineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 107bbd5d0b3SPeter Brune if (changed_y) { 108bbd5d0b3SPeter Brune ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 109bbd5d0b3SPeter Brune } else { 110bbd5d0b3SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 111bbd5d0b3SPeter Brune } 112bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 113*6a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 114*6a388c36SPeter Brune if (domainerror) { 115*6a388c36SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 116bbd5d0b3SPeter Brune PetscFunctionReturn(0); 117bbd5d0b3SPeter Brune } 118bbd5d0b3SPeter Brune 119*6a388c36SPeter Brune ierr = LineSearchComputeNorms(linesearch);CHKERRQ(ierr); 120*6a388c36SPeter Brune ierr = LineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 121bbd5d0b3SPeter Brune 122*6a388c36SPeter Brune ierr = LineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 123*6a388c36SPeter Brune 124*6a388c36SPeter Brune if (monitor) { 125*6a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 126*6a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", lambda, gnorm);CHKERRQ(ierr); 127*6a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 128bbd5d0b3SPeter Brune } 129*6a388c36SPeter Brune if (lambda <= steptol) { 130*6a388c36SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 131bbd5d0b3SPeter Brune } 132bbd5d0b3SPeter Brune PetscFunctionReturn(0); 133bbd5d0b3SPeter Brune } 134bbd5d0b3SPeter Brune 135bbd5d0b3SPeter Brune EXTERN_C_BEGIN 136bbd5d0b3SPeter Brune #undef __FUNCT__ 137bbd5d0b3SPeter Brune #define __FUNCT__ "LineSearchCreate_CP" 138bbd5d0b3SPeter Brune PetscErrorCode LineSearchCreate_CP(LineSearch linesearch) 139bbd5d0b3SPeter Brune { 140bbd5d0b3SPeter Brune PetscFunctionBegin; 141bbd5d0b3SPeter Brune linesearch->ops->apply = LineSearchApply_CP; 142bbd5d0b3SPeter Brune linesearch->ops->destroy = PETSC_NULL; 143bbd5d0b3SPeter Brune linesearch->ops->setfromoptions = PETSC_NULL; 144bbd5d0b3SPeter Brune linesearch->ops->reset = PETSC_NULL; 145bbd5d0b3SPeter Brune linesearch->ops->view = PETSC_NULL; 146bbd5d0b3SPeter Brune linesearch->ops->setup = PETSC_NULL; 147bbd5d0b3SPeter Brune PetscFunctionReturn(0); 148bbd5d0b3SPeter Brune } 149bbd5d0b3SPeter Brune EXTERN_C_END 150