1bbd5d0b3SPeter Brune #include <private/linesearchimpl.h> 26a388c36SPeter Brune #include <petscsnes.h> 3bbd5d0b3SPeter Brune 4bbd5d0b3SPeter Brune #undef __FUNCT__ 5b9423516SPeter Brune #define __FUNCT__ "PetscLineSearchApply_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 @*/ 406188f407SPeter Brune PetscErrorCode PetscLineSearchApply_CP(PetscLineSearch linesearch) 41bbd5d0b3SPeter Brune { 42bbd5d0b3SPeter Brune 436a388c36SPeter Brune PetscBool changed_y, changed_w, domainerror; 44bbd5d0b3SPeter Brune PetscErrorCode ierr; 456a388c36SPeter Brune Vec X, Y, F, W; 466a388c36SPeter Brune SNES snes; 47516fe3c3SPeter Brune PetscReal xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep; 48bbd5d0b3SPeter Brune 49bbd5d0b3SPeter Brune PetscReal lambda, lambda_old, lambda_update, delLambda; 50bbd5d0b3SPeter Brune PetscScalar fty, fty_old; 516a388c36SPeter Brune PetscInt i, max_its; 526a388c36SPeter Brune 536a388c36SPeter Brune PetscViewer monitor; 54bbd5d0b3SPeter Brune 55bbd5d0b3SPeter Brune PetscFunctionBegin; 566a388c36SPeter Brune 576188f407SPeter Brune ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr); 586188f407SPeter Brune ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 596188f407SPeter Brune ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 606188f407SPeter Brune ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 616188f407SPeter Brune ierr = PetscLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its);CHKERRQ(ierr); 626188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 636188f407SPeter Brune ierr = PetscLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 646a388c36SPeter Brune 65bbd5d0b3SPeter Brune /* precheck */ 666188f407SPeter Brune ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 67bbd5d0b3SPeter Brune lambda_old = 0.0; 68bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr); 696a388c36SPeter Brune for (i = 0; i < max_its; i++) { 70bbd5d0b3SPeter Brune 71bbd5d0b3SPeter Brune /* compute the norm at lambda */ 72bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 73bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 74*9bd66eb0SPeter Brune if (linesearch->ops->viproject) { 75*9bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 76*9bd66eb0SPeter Brune } 77bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 78bbd5d0b3SPeter Brune 79bbd5d0b3SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 80bbd5d0b3SPeter Brune 81bbd5d0b3SPeter Brune delLambda = lambda - lambda_old; 826a388c36SPeter Brune if (PetscAbsReal(delLambda) < steptol) break; 836a388c36SPeter Brune if (monitor) { 846a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 856a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g], ftys = [%g, %g]\n", 86bbd5d0b3SPeter Brune lambda, lambda_old, PetscRealPart(fty), PetscRealPart(fty_old));CHKERRQ(ierr); 876a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 88bbd5d0b3SPeter Brune } 89bbd5d0b3SPeter Brune 90bbd5d0b3SPeter Brune /* compute the search direction */ 91bbd5d0b3SPeter Brune lambda_update = PetscRealPart((fty*lambda_old - fty_old*lambda) / (fty - fty_old)); 92bbd5d0b3SPeter Brune if (PetscIsInfOrNanScalar(lambda_update)) break; 936a388c36SPeter Brune if (lambda_update > maxstep) { 94bbd5d0b3SPeter Brune break; 95bbd5d0b3SPeter Brune } 96bbd5d0b3SPeter Brune 97bbd5d0b3SPeter Brune /* compute the new state of the line search */ 98bbd5d0b3SPeter Brune lambda_old = lambda; 99bbd5d0b3SPeter Brune lambda = lambda_update; 100bbd5d0b3SPeter Brune fty_old = fty; 101bbd5d0b3SPeter Brune } 102bbd5d0b3SPeter Brune /* construct the solution */ 103bbd5d0b3SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 104bbd5d0b3SPeter Brune ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 105*9bd66eb0SPeter Brune if (linesearch->ops->viproject) { 106*9bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 107*9bd66eb0SPeter Brune } 108bbd5d0b3SPeter Brune /* postcheck */ 1096188f407SPeter Brune ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 110bbd5d0b3SPeter Brune if (changed_y) { 111bbd5d0b3SPeter Brune ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 112*9bd66eb0SPeter Brune if (linesearch->ops->viproject) { 113*9bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); 114*9bd66eb0SPeter Brune } 115bbd5d0b3SPeter Brune } else { 116bbd5d0b3SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 117bbd5d0b3SPeter Brune } 118bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1196a388c36SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1206a388c36SPeter Brune if (domainerror) { 1216188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 122bbd5d0b3SPeter Brune PetscFunctionReturn(0); 123bbd5d0b3SPeter Brune } 124bbd5d0b3SPeter Brune 1256188f407SPeter Brune ierr = PetscLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 1266188f407SPeter Brune ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 127bbd5d0b3SPeter Brune 1286188f407SPeter Brune ierr = PetscLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 1296a388c36SPeter Brune 1306a388c36SPeter Brune if (monitor) { 1316a388c36SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 1326a388c36SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", lambda, gnorm);CHKERRQ(ierr); 1336a388c36SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 134bbd5d0b3SPeter Brune } 1356a388c36SPeter Brune if (lambda <= steptol) { 1366188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 137bbd5d0b3SPeter Brune } 138bbd5d0b3SPeter Brune PetscFunctionReturn(0); 139bbd5d0b3SPeter Brune } 140bbd5d0b3SPeter Brune 141bbd5d0b3SPeter Brune EXTERN_C_BEGIN 142bbd5d0b3SPeter Brune #undef __FUNCT__ 1436188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate_CP" 1446188f407SPeter Brune PetscErrorCode PetscLineSearchCreate_CP(PetscLineSearch linesearch) 145bbd5d0b3SPeter Brune { 146bbd5d0b3SPeter Brune PetscFunctionBegin; 1476188f407SPeter Brune linesearch->ops->apply = PetscLineSearchApply_CP; 148bbd5d0b3SPeter Brune linesearch->ops->destroy = PETSC_NULL; 149bbd5d0b3SPeter Brune linesearch->ops->setfromoptions = PETSC_NULL; 150bbd5d0b3SPeter Brune linesearch->ops->reset = PETSC_NULL; 151bbd5d0b3SPeter Brune linesearch->ops->view = PETSC_NULL; 152bbd5d0b3SPeter Brune linesearch->ops->setup = PETSC_NULL; 153bbd5d0b3SPeter Brune PetscFunctionReturn(0); 154bbd5d0b3SPeter Brune } 155bbd5d0b3SPeter Brune EXTERN_C_END 156