1b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2bf7f4e0aSPeter Brune 3f1c6b773SPeter Brune PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 4f1c6b773SPeter Brune PetscFList SNESLineSearchList = PETSC_NULL; 5bf7f4e0aSPeter Brune 6f1c6b773SPeter Brune PetscClassId SNESLINESEARCH_CLASSID; 7f1c6b773SPeter Brune PetscLogEvent SNESLineSearch_Apply; 8bf7f4e0aSPeter Brune 9bf7f4e0aSPeter Brune #undef __FUNCT__ 10f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate" 11f40b411bSPeter Brune /*@ 12cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 13f40b411bSPeter Brune 14cd7522eaSPeter Brune Logically Collective on Comm 15f40b411bSPeter Brune 16f40b411bSPeter Brune Input Parameters: 17cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 18f40b411bSPeter Brune 19f40b411bSPeter Brune Output Parameters: 20f40b411bSPeter Brune . outlinesearch - the line search instance. 21f40b411bSPeter Brune 22*78bcb3b5SPeter Brune Level: beginner 23f40b411bSPeter Brune 24cd7522eaSPeter Brune .keywords: LineSearch, create, context 25f40b411bSPeter Brune 26f40b411bSPeter Brune .seealso: LineSearchDestroy() 27f40b411bSPeter Brune @*/ 28f40b411bSPeter Brune 29f1c6b773SPeter Brune PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) { 30bf7f4e0aSPeter Brune PetscErrorCode ierr; 31f1c6b773SPeter Brune SNESLineSearch linesearch; 32bf7f4e0aSPeter Brune PetscFunctionBegin; 33ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 34ea5d4fccSPeter Brune *outlinesearch = PETSC_NULL; 35f1c6b773SPeter Brune ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, 0, 36f1c6b773SPeter Brune "SNESLineSearch","Line-search method","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 37bf7f4e0aSPeter Brune 38bf7f4e0aSPeter Brune linesearch->ops->precheckstep = PETSC_NULL; 39bf7f4e0aSPeter Brune linesearch->ops->postcheckstep = PETSC_NULL; 40bf7f4e0aSPeter Brune 419bd66eb0SPeter Brune linesearch->vec_sol_new = PETSC_NULL; 429bd66eb0SPeter Brune linesearch->vec_func_new = PETSC_NULL; 439bd66eb0SPeter Brune linesearch->vec_sol = PETSC_NULL; 449bd66eb0SPeter Brune linesearch->vec_func = PETSC_NULL; 459bd66eb0SPeter Brune linesearch->vec_update = PETSC_NULL; 469bd66eb0SPeter Brune 47bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 48bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 49bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 50bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 51bf7f4e0aSPeter Brune linesearch->success = PETSC_TRUE; 52bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 53bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 54bf7f4e0aSPeter Brune linesearch->damping = 1.0; 55bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 56bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 57516fe3c3SPeter Brune linesearch->rtol = 1e-8; 58516fe3c3SPeter Brune linesearch->atol = 1e-15; 59516fe3c3SPeter Brune linesearch->ltol = 1e-8; 60bf7f4e0aSPeter Brune linesearch->precheckctx = PETSC_NULL; 61bf7f4e0aSPeter Brune linesearch->postcheckctx = PETSC_NULL; 6259405d5eSPeter Brune linesearch->max_its = 1; 63bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 64bf7f4e0aSPeter Brune *outlinesearch = linesearch; 65bf7f4e0aSPeter Brune PetscFunctionReturn(0); 66bf7f4e0aSPeter Brune } 67bf7f4e0aSPeter Brune 68bf7f4e0aSPeter Brune #undef __FUNCT__ 69f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp" 70f40b411bSPeter Brune /*@ 71*78bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 72*78bcb3b5SPeter Brune any required vectors. 73f40b411bSPeter Brune 74cd7522eaSPeter Brune Collective on SNESLineSearch 75f40b411bSPeter Brune 76f40b411bSPeter Brune Input Parameters: 77f40b411bSPeter Brune . linesearch - The LineSearch instance. 78f40b411bSPeter Brune 79cd7522eaSPeter Brune Notes: 80cd7522eaSPeter Brune For most cases, this needn't be called outside of SNESLineSearchApply(). 81cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 82*78bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 83cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 84cd7522eaSPeter Brune allocated upfront. 85cd7522eaSPeter Brune 86cd7522eaSPeter Brune 87*78bcb3b5SPeter Brune Level: advanced 88f40b411bSPeter Brune 89f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 90f40b411bSPeter Brune 91f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 92f40b411bSPeter Brune @*/ 93f40b411bSPeter Brune 94f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) { 95bf7f4e0aSPeter Brune PetscErrorCode ierr; 96bf7f4e0aSPeter Brune PetscFunctionBegin; 97bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 98f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNES_LINESEARCH_BASIC);CHKERRQ(ierr); 99bf7f4e0aSPeter Brune } 100bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 1019bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 102bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 1039bd66eb0SPeter Brune } 1049bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 1059bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 1069bd66eb0SPeter Brune } 107bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 108bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 109bf7f4e0aSPeter Brune } 110bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 111bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 112bf7f4e0aSPeter Brune } 113bf7f4e0aSPeter Brune PetscFunctionReturn(0); 114bf7f4e0aSPeter Brune } 115bf7f4e0aSPeter Brune 116bf7f4e0aSPeter Brune #undef __FUNCT__ 117f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset" 118f40b411bSPeter Brune 119f40b411bSPeter Brune /*@ 120*78bcb3b5SPeter Brune SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search. 121f40b411bSPeter Brune 122f1c6b773SPeter Brune Collective on SNESLineSearch 123f40b411bSPeter Brune 124f40b411bSPeter Brune Input Parameters: 125f40b411bSPeter Brune . linesearch - The LineSearch instance. 126f40b411bSPeter Brune 127*78bcb3b5SPeter Brune Level: intermediate 128f40b411bSPeter Brune 129cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 130f40b411bSPeter Brune 131f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 132f40b411bSPeter Brune @*/ 133f40b411bSPeter Brune 134f1c6b773SPeter Brune PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) { 135bf7f4e0aSPeter Brune PetscErrorCode ierr; 136bf7f4e0aSPeter Brune PetscFunctionBegin; 137bf7f4e0aSPeter Brune if (linesearch->ops->reset) { 138bf7f4e0aSPeter Brune (*linesearch->ops->reset)(linesearch); 139bf7f4e0aSPeter Brune } 140bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 141bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 142bf7f4e0aSPeter Brune 143bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 144bf7f4e0aSPeter Brune linesearch->nwork = 0; 145bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 146bf7f4e0aSPeter Brune PetscFunctionReturn(0); 147bf7f4e0aSPeter Brune } 148bf7f4e0aSPeter Brune 14986d74e61SPeter Brune 15086d74e61SPeter Brune #undef __FUNCT__ 151f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck" 15286d74e61SPeter Brune /*@C 153f1c6b773SPeter Brune SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine. 15486d74e61SPeter Brune 155f1c6b773SPeter Brune Logically Collective on SNESLineSearch 15686d74e61SPeter Brune 15786d74e61SPeter Brune Input Parameters: 158f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 15986d74e61SPeter Brune . func - [optional] function evaluation routine 16086d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 16186d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 16286d74e61SPeter Brune 16386d74e61SPeter Brune Calling sequence of func: 164f1c6b773SPeter Brune $ func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 16586d74e61SPeter Brune 16686d74e61SPeter Brune + x - solution vector 16786d74e61SPeter Brune . y - search direction vector 168cd7522eaSPeter Brune - changed - flag to indicate the precheck changed x or y. 16986d74e61SPeter Brune 17086d74e61SPeter Brune Level: intermediate 17186d74e61SPeter Brune 17286d74e61SPeter Brune .keywords: set, linesearch, pre-check 17386d74e61SPeter Brune 174f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck() 17586d74e61SPeter Brune @*/ 176f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx) 17786d74e61SPeter Brune { 1789bd66eb0SPeter Brune PetscFunctionBegin; 179f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 18086d74e61SPeter Brune if (func) linesearch->ops->precheckstep = func; 18186d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 18286d74e61SPeter Brune PetscFunctionReturn(0); 18386d74e61SPeter Brune } 18486d74e61SPeter Brune 18586d74e61SPeter Brune 18686d74e61SPeter Brune #undef __FUNCT__ 187f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck" 18886d74e61SPeter Brune /*@C 189cd7522eaSPeter Brune SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine. 19086d74e61SPeter Brune 19186d74e61SPeter Brune Input Parameters: 192f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 19386d74e61SPeter Brune 19486d74e61SPeter Brune Output Parameters: 19586d74e61SPeter Brune + func - [optional] function evaluation routine 19686d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 19786d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 19886d74e61SPeter Brune 19986d74e61SPeter Brune Level: intermediate 20086d74e61SPeter Brune 20186d74e61SPeter Brune .keywords: get, linesearch, pre-check 20286d74e61SPeter Brune 203f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 20486d74e61SPeter Brune @*/ 205f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx) 20686d74e61SPeter Brune { 2079bd66eb0SPeter Brune PetscFunctionBegin; 208f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 20986d74e61SPeter Brune if (func) *func = linesearch->ops->precheckstep; 21086d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 21186d74e61SPeter Brune PetscFunctionReturn(0); 21286d74e61SPeter Brune } 21386d74e61SPeter Brune 21486d74e61SPeter Brune 21586d74e61SPeter Brune #undef __FUNCT__ 216f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck" 21786d74e61SPeter Brune /*@C 218f1c6b773SPeter Brune SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine. 21986d74e61SPeter Brune 220f1c6b773SPeter Brune Logically Collective on SNESLineSearch 22186d74e61SPeter Brune 22286d74e61SPeter Brune Input Parameters: 223f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 22486d74e61SPeter Brune . func - [optional] function evaluation routine 22586d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 22686d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 22786d74e61SPeter Brune 22886d74e61SPeter Brune Calling sequence of func: 229f1c6b773SPeter Brune $ func (SNESLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y); 23086d74e61SPeter Brune 23186d74e61SPeter Brune + x - old solution vector 23286d74e61SPeter Brune . y - search direction vector 23386d74e61SPeter Brune . w - new solution vector 23486d74e61SPeter Brune . changed_y - indicates that the line search changed y. 23586d74e61SPeter Brune . changed_w - indicates that the line search changed w. 23686d74e61SPeter Brune 23786d74e61SPeter Brune Level: intermediate 23886d74e61SPeter Brune 23986d74e61SPeter Brune .keywords: set, linesearch, post-check 24086d74e61SPeter Brune 241f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 24286d74e61SPeter Brune @*/ 243f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx) 24486d74e61SPeter Brune { 24586d74e61SPeter Brune PetscFunctionBegin; 246f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 24786d74e61SPeter Brune if (func) linesearch->ops->postcheckstep = func; 24886d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 24986d74e61SPeter Brune PetscFunctionReturn(0); 25086d74e61SPeter Brune } 25186d74e61SPeter Brune 25286d74e61SPeter Brune 25386d74e61SPeter Brune #undef __FUNCT__ 254f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck" 25586d74e61SPeter Brune /*@C 256f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 25786d74e61SPeter Brune 25886d74e61SPeter Brune Input Parameters: 259f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 26086d74e61SPeter Brune 26186d74e61SPeter Brune Output Parameters: 26286d74e61SPeter Brune + func - [optional] function evaluation routine 26386d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 26486d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 26586d74e61SPeter Brune 26686d74e61SPeter Brune Level: intermediate 26786d74e61SPeter Brune 26886d74e61SPeter Brune .keywords: get, linesearch, post-check 26986d74e61SPeter Brune 270f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck() 27186d74e61SPeter Brune @*/ 272f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc *func,void **ctx) 27386d74e61SPeter Brune { 2749bd66eb0SPeter Brune PetscFunctionBegin; 275f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 27686d74e61SPeter Brune if (func) *func = linesearch->ops->postcheckstep; 27786d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 27886d74e61SPeter Brune PetscFunctionReturn(0); 27986d74e61SPeter Brune } 28086d74e61SPeter Brune 28186d74e61SPeter Brune 282bf7f4e0aSPeter Brune #undef __FUNCT__ 283f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck" 284f40b411bSPeter Brune /*@ 285f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 286f40b411bSPeter Brune 287cd7522eaSPeter Brune Logically Collective on SNESLineSearch 288f40b411bSPeter Brune 289f40b411bSPeter Brune Input Parameters: 290f40b411bSPeter Brune . linesearch - The linesearch instance. 291f40b411bSPeter Brune 292f40b411bSPeter Brune Output Parameters: 293*78bcb3b5SPeter Brune . changed - Indicator that the precheck routine has changed anything. 294f40b411bSPeter Brune 295f40b411bSPeter Brune Level: Beginner 296f40b411bSPeter Brune 297f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 298f40b411bSPeter Brune 299f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 300f40b411bSPeter Brune @*/ 301f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, PetscBool * changed) 302bf7f4e0aSPeter Brune { 303bf7f4e0aSPeter Brune PetscErrorCode ierr; 304bf7f4e0aSPeter Brune PetscFunctionBegin; 305bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 306bf7f4e0aSPeter Brune if (linesearch->ops->precheckstep) { 307c87759e9SPeter Brune ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed, linesearch->precheckctx);CHKERRQ(ierr); 308bf7f4e0aSPeter Brune } 309bf7f4e0aSPeter Brune PetscFunctionReturn(0); 310bf7f4e0aSPeter Brune } 311bf7f4e0aSPeter Brune 312bf7f4e0aSPeter Brune #undef __FUNCT__ 313f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck" 314f40b411bSPeter Brune /*@ 315f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 316f40b411bSPeter Brune 317cd7522eaSPeter Brune Logically Collective on SNESLineSearch 318f40b411bSPeter Brune 319f40b411bSPeter Brune Input Parameters: 320f40b411bSPeter Brune . linesearch - The linesearch instance. 321f40b411bSPeter Brune 322f40b411bSPeter Brune Output Parameters: 323*78bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 324*78bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 325f40b411bSPeter Brune 326f40b411bSPeter Brune Level: Intermediate 327f40b411bSPeter Brune 328f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 329f40b411bSPeter Brune 330f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 331f40b411bSPeter Brune @*/ 332f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, PetscBool * changed_Y, PetscBool * changed_W) 333bf7f4e0aSPeter Brune { 334bf7f4e0aSPeter Brune PetscErrorCode ierr; 335bf7f4e0aSPeter Brune PetscFunctionBegin; 336bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 337bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 338bf7f4e0aSPeter Brune if (linesearch->ops->postcheckstep) { 339c87759e9SPeter Brune ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, linesearch->vec_sol_new, changed_Y, changed_W, linesearch->postcheckctx);CHKERRQ(ierr); 34086d74e61SPeter Brune } 34186d74e61SPeter Brune PetscFunctionReturn(0); 34286d74e61SPeter Brune } 34386d74e61SPeter Brune 34486d74e61SPeter Brune 34586d74e61SPeter Brune #undef __FUNCT__ 346f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard" 34786d74e61SPeter Brune /*@C 34886d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 34986d74e61SPeter Brune 350cd7522eaSPeter Brune Logically Collective on SNESLineSearch 35186d74e61SPeter Brune 35286d74e61SPeter Brune Input Arguments: 35386d74e61SPeter Brune + linesearch - linesearch context 35486d74e61SPeter Brune . X - base state for this step 35586d74e61SPeter Brune . Y - initial correction 35686d74e61SPeter Brune 35786d74e61SPeter Brune Output Arguments: 35886d74e61SPeter Brune + Y - correction, possibly modified 35986d74e61SPeter Brune - changed - flag indicating that Y was modified 36086d74e61SPeter Brune 36186d74e61SPeter Brune Options Database Key: 362cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 363cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 36486d74e61SPeter Brune 36586d74e61SPeter Brune Level: advanced 36686d74e61SPeter Brune 36786d74e61SPeter Brune Notes: 36886d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 36986d74e61SPeter Brune 37086d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 37186d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 37286d74e61SPeter Brune is generally not useful when using a Newton linearization. 37386d74e61SPeter Brune 37486d74e61SPeter Brune Reference: 37586d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 37686d74e61SPeter Brune 37786d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 37886d74e61SPeter Brune @*/ 379f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 38086d74e61SPeter Brune { 38186d74e61SPeter Brune PetscErrorCode ierr; 38286d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 38386d74e61SPeter Brune Vec Ylast; 38486d74e61SPeter Brune PetscScalar dot; 385c87759e9SPeter Brune 38686d74e61SPeter Brune PetscInt iter; 38786d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 38886d74e61SPeter Brune SNES snes; 38986d74e61SPeter Brune 39086d74e61SPeter Brune PetscFunctionBegin; 391f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 39286d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 39386d74e61SPeter Brune if (!Ylast) { 39486d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 39586d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 39686d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 39786d74e61SPeter Brune } 39886d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 39986d74e61SPeter Brune if (iter < 2) { 40086d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 40186d74e61SPeter Brune *changed = PETSC_FALSE; 40286d74e61SPeter Brune PetscFunctionReturn(0); 40386d74e61SPeter Brune } 40486d74e61SPeter Brune 40586d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 40686d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 40786d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 40886d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 40986d74e61SPeter Brune theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 41086d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 41186d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 41286d74e61SPeter Brune /* Modify the step Y */ 41386d74e61SPeter Brune PetscReal alpha,ydiffnorm; 41486d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 41586d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 41686d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 41786d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 41886d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 41986d74e61SPeter Brune ierr = PetscInfo3(snes,"Angle %G degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,angle,alpha);CHKERRQ(ierr); 42086d74e61SPeter Brune } else { 42186d74e61SPeter Brune ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr); 42286d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 42386d74e61SPeter Brune *changed = PETSC_FALSE; 424bf7f4e0aSPeter Brune } 425bf7f4e0aSPeter Brune PetscFunctionReturn(0); 426bf7f4e0aSPeter Brune } 427bf7f4e0aSPeter Brune 428bf7f4e0aSPeter Brune #undef __FUNCT__ 429f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply" 430f40b411bSPeter Brune /*@ 431cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 432f40b411bSPeter Brune 433f1c6b773SPeter Brune Collective on SNESLineSearch 434f40b411bSPeter Brune 435f40b411bSPeter Brune Input Parameters: 436f40b411bSPeter Brune + linesearch - The linesearch instance. 437f40b411bSPeter Brune . X - The current solution. 438f40b411bSPeter Brune . F - The current function. 439f40b411bSPeter Brune . fnorm - The current norm. 440f40b411bSPeter Brune . Y - The search direction. 441f40b411bSPeter Brune 442f40b411bSPeter Brune Output Parameters: 443f40b411bSPeter Brune + X - The new solution. 444f40b411bSPeter Brune . F - The new function. 445f40b411bSPeter Brune - fnorm - The new function norm. 446f40b411bSPeter Brune 447cd7522eaSPeter Brune Options Database Keys: 448cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method 449cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 450cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter. 451cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 452cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 453cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches. 454cd7522eaSPeter Brune 455cd7522eaSPeter Brune Notes: 456cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 457cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 458cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 459cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 460cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 461cd7522eaSPeter Brune therefore may be fairly expensive. 462cd7522eaSPeter Brune 463f40b411bSPeter Brune Level: Intermediate 464f40b411bSPeter Brune 465f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 466f40b411bSPeter Brune 467cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction() 468f40b411bSPeter Brune @*/ 469f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) { 470bf7f4e0aSPeter Brune PetscErrorCode ierr; 471bf7f4e0aSPeter Brune PetscFunctionBegin; 472bf7f4e0aSPeter Brune 473bf7f4e0aSPeter Brune /* check the pointers */ 474f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 475bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 476bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 477bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 478bf7f4e0aSPeter Brune 479bf7f4e0aSPeter Brune linesearch->success = PETSC_TRUE; 480bf7f4e0aSPeter Brune 481bf7f4e0aSPeter Brune linesearch->vec_sol = X; 482bf7f4e0aSPeter Brune linesearch->vec_update = Y; 483bf7f4e0aSPeter Brune linesearch->vec_func = F; 484bf7f4e0aSPeter Brune 485f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 486bf7f4e0aSPeter Brune 487bf7f4e0aSPeter Brune if (!linesearch->keeplambda) 488bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 489bf7f4e0aSPeter Brune 490bf7f4e0aSPeter Brune if (fnorm) { 491bf7f4e0aSPeter Brune linesearch->fnorm = *fnorm; 492bf7f4e0aSPeter Brune } else { 493bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 494bf7f4e0aSPeter Brune } 495bf7f4e0aSPeter Brune 496f1c6b773SPeter Brune ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 497bf7f4e0aSPeter Brune 498bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 499bf7f4e0aSPeter Brune 500f1c6b773SPeter Brune ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 501bf7f4e0aSPeter Brune 502bf7f4e0aSPeter Brune if (fnorm) 503bf7f4e0aSPeter Brune *fnorm = linesearch->fnorm; 504bf7f4e0aSPeter Brune PetscFunctionReturn(0); 505bf7f4e0aSPeter Brune } 506bf7f4e0aSPeter Brune 507bf7f4e0aSPeter Brune #undef __FUNCT__ 508f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy" 509f40b411bSPeter Brune /*@ 510f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 511f40b411bSPeter Brune 512f1c6b773SPeter Brune Collective on SNESLineSearch 513f40b411bSPeter Brune 514f40b411bSPeter Brune Input Parameters: 515f40b411bSPeter Brune . linesearch - The linesearch instance. 516f40b411bSPeter Brune 517f40b411bSPeter Brune Level: Intermediate 518f40b411bSPeter Brune 519*78bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 520f40b411bSPeter Brune 521cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 522f40b411bSPeter Brune @*/ 523f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) { 524bf7f4e0aSPeter Brune PetscErrorCode ierr; 525bf7f4e0aSPeter Brune PetscFunctionBegin; 526bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 527f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 528bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 529bf7f4e0aSPeter Brune ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr); 530f1c6b773SPeter Brune ierr = SNESLineSearchReset(*linesearch); 531bf7f4e0aSPeter Brune if ((*linesearch)->ops->destroy) { 532bf7f4e0aSPeter Brune (*linesearch)->ops->destroy(*linesearch); 533bf7f4e0aSPeter Brune } 534bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 535e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 536bf7f4e0aSPeter Brune PetscFunctionReturn(0); 537bf7f4e0aSPeter Brune } 538bf7f4e0aSPeter Brune 539bf7f4e0aSPeter Brune #undef __FUNCT__ 540f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor" 541f40b411bSPeter Brune /*@ 542cd7522eaSPeter Brune SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search. 543bf7f4e0aSPeter Brune 544bf7f4e0aSPeter Brune Input Parameters: 545bf7f4e0aSPeter Brune + snes - nonlinear context obtained from SNESCreate() 546bf7f4e0aSPeter Brune - flg - PETSC_TRUE to monitor the line search 547bf7f4e0aSPeter Brune 548bf7f4e0aSPeter Brune Logically Collective on SNES 549bf7f4e0aSPeter Brune 550bf7f4e0aSPeter Brune Options Database: 551cd7522eaSPeter Brune . -snes_linesearch_monitor - enables the monitor. 552bf7f4e0aSPeter Brune 553bf7f4e0aSPeter Brune Level: intermediate 554bf7f4e0aSPeter Brune 555bf7f4e0aSPeter Brune 556cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer 557bf7f4e0aSPeter Brune @*/ 558f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg) 559bf7f4e0aSPeter Brune { 560bf7f4e0aSPeter Brune 561bf7f4e0aSPeter Brune PetscErrorCode ierr; 562bf7f4e0aSPeter Brune PetscFunctionBegin; 563bf7f4e0aSPeter Brune if (flg && !linesearch->monitor) { 564bf7f4e0aSPeter Brune ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr); 565bf7f4e0aSPeter Brune } else if (!flg && linesearch->monitor) { 566bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 567bf7f4e0aSPeter Brune } 568bf7f4e0aSPeter Brune PetscFunctionReturn(0); 569bf7f4e0aSPeter Brune } 570bf7f4e0aSPeter Brune 571bf7f4e0aSPeter Brune #undef __FUNCT__ 572f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor" 573f40b411bSPeter Brune /*@ 574cd7522eaSPeter Brune SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor. 5756a388c36SPeter Brune 576f40b411bSPeter Brune Input Parameters: 577f40b411bSPeter Brune . linesearch - linesearch context. 578f40b411bSPeter Brune 579f40b411bSPeter Brune Input Parameters: 580f40b411bSPeter Brune . monitor - monitor context. 581f40b411bSPeter Brune 582f40b411bSPeter Brune Logically Collective on SNES 583f40b411bSPeter Brune 584f40b411bSPeter Brune 585f40b411bSPeter Brune Options Database Keys: 586cd7522eaSPeter Brune . -snes_linesearch_monitor - enables the monitor. 587f40b411bSPeter Brune 588f40b411bSPeter Brune Level: intermediate 589f40b411bSPeter Brune 590f40b411bSPeter Brune 591cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer 592f40b411bSPeter Brune @*/ 593f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 5946a388c36SPeter Brune { 5956a388c36SPeter Brune 5966a388c36SPeter Brune PetscFunctionBegin; 597f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 5986a388c36SPeter Brune if (monitor) { 5996a388c36SPeter Brune PetscValidPointer(monitor, 2); 6006a388c36SPeter Brune *monitor = linesearch->monitor; 6016a388c36SPeter Brune } 6026a388c36SPeter Brune PetscFunctionReturn(0); 6036a388c36SPeter Brune } 6046a388c36SPeter Brune 6056a388c36SPeter Brune #undef __FUNCT__ 606f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions" 607f40b411bSPeter Brune /*@ 608f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 609f40b411bSPeter Brune 610f40b411bSPeter Brune Input Parameters: 611f40b411bSPeter Brune . linesearch - linesearch context. 612f40b411bSPeter Brune 613f40b411bSPeter Brune Options Database Keys: 614cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method 615cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 616cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter. 617cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 618cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 61959405d5eSPeter Brune . -snes_linesearch_order - The polynomial order of approximation used in the linesearch 620cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches. 621f40b411bSPeter Brune 622f1c6b773SPeter Brune Logically Collective on SNESLineSearch 623f40b411bSPeter Brune 624f40b411bSPeter Brune Level: intermediate 625f40b411bSPeter Brune 626f40b411bSPeter Brune 627f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 628f40b411bSPeter Brune @*/ 629f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) { 630bf7f4e0aSPeter Brune PetscErrorCode ierr; 63159405d5eSPeter Brune const char *orders[] = {"linear", "quadratic", "cubic"}; 63259405d5eSPeter Brune PetscInt indx = 0; 633f1c6b773SPeter Brune const char *deft = SNES_LINESEARCH_BASIC; 634bf7f4e0aSPeter Brune char type[256]; 635bf7f4e0aSPeter Brune PetscBool flg, set; 636bf7f4e0aSPeter Brune PetscFunctionBegin; 637f1c6b773SPeter Brune if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 638bf7f4e0aSPeter Brune 639bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 640bf7f4e0aSPeter Brune if (((PetscObject)linesearch)->type_name) { 641bf7f4e0aSPeter Brune deft = ((PetscObject)linesearch)->type_name; 642bf7f4e0aSPeter Brune } 6437a35526eSPeter Brune ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 644bf7f4e0aSPeter Brune if (flg) { 645f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 646bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 647f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 648bf7f4e0aSPeter Brune } 649bf7f4e0aSPeter Brune if (linesearch->ops->setfromoptions) { 650bf7f4e0aSPeter Brune (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr); 651bf7f4e0aSPeter Brune } 652bf7f4e0aSPeter Brune 6537a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor", 654bf7f4e0aSPeter Brune linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 655f1c6b773SPeter Brune if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 656bf7f4e0aSPeter Brune 6577a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr); 6587a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_rtol","Tolerance for iterative line search","SNESLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr); 6597a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr); 6607a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr); 6617a35526eSPeter Brune ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr); 6627a35526eSPeter Brune 6637a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 6647a35526eSPeter Brune if (set) { 6657a35526eSPeter Brune if (flg) { 6667a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 6677a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 6687a35526eSPeter Brune "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr); 6697a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 6707a35526eSPeter Brune } else { 6717a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 6727a35526eSPeter Brune } 6737a35526eSPeter Brune } 67459405d5eSPeter Brune ierr = PetscOptionsEList("-snes_linesearch_order", "Order of approximation", "SNESLineSearch", orders,3,orders[(PetscInt)linesearch->order],&indx,&flg);CHKERRQ(ierr); 67559405d5eSPeter Brune if (flg) { 67659405d5eSPeter Brune switch (indx) { 67759405d5eSPeter Brune case 0: linesearch->order = SNES_LINESEARCH_LINEAR; 67859405d5eSPeter Brune break; 67959405d5eSPeter Brune case 1: linesearch->order = SNES_LINESEARCH_QUADRATIC; 68059405d5eSPeter Brune break; 68159405d5eSPeter Brune case 2: linesearch->order = SNES_LINESEARCH_CUBIC; 68259405d5eSPeter Brune break; 68359405d5eSPeter Brune } 68459405d5eSPeter Brune } 68559405d5eSPeter Brune 6867a35526eSPeter Brune 687bf7f4e0aSPeter Brune ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr); 688bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 689bf7f4e0aSPeter Brune PetscFunctionReturn(0); 690bf7f4e0aSPeter Brune } 691bf7f4e0aSPeter Brune 692bf7f4e0aSPeter Brune #undef __FUNCT__ 693f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView" 694f40b411bSPeter Brune /*@ 695cd7522eaSPeter Brune SNESLineSearchView - Prints useful information about the line search not 696cd7522eaSPeter Brune related to an individual call. 697f40b411bSPeter Brune 698f40b411bSPeter Brune Input Parameters: 699f40b411bSPeter Brune . linesearch - linesearch context. 700f40b411bSPeter Brune 701f1c6b773SPeter Brune Logically Collective on SNESLineSearch 702f40b411bSPeter Brune 703f40b411bSPeter Brune Level: intermediate 704f40b411bSPeter Brune 705f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 706f40b411bSPeter Brune @*/ 7077f1410a3SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) { 7087f1410a3SPeter Brune PetscErrorCode ierr; 7097f1410a3SPeter Brune PetscBool iascii; 710bf7f4e0aSPeter Brune PetscFunctionBegin; 7117f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 7127f1410a3SPeter Brune if (!viewer) { 7137f1410a3SPeter Brune ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr); 7147f1410a3SPeter Brune } 7157f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 7167f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 717f40b411bSPeter Brune 7187f1410a3SPeter Brune ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7197f1410a3SPeter Brune if (iascii) { 7207f1410a3SPeter Brune ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr); 7217f1410a3SPeter Brune if (linesearch->ops->view) { 7227f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7237f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 7247f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7257f1410a3SPeter Brune } 7267f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maxstep=%G, minlambda=%G\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr); 7277f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, lambda=%G\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr); 7287f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 7297f1410a3SPeter Brune if (linesearch->ops->precheckstep) { 7307f1410a3SPeter Brune if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) { 7317f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 7327f1410a3SPeter Brune } else { 7337f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 7347f1410a3SPeter Brune } 7357f1410a3SPeter Brune } 7367f1410a3SPeter Brune if (linesearch->ops->postcheckstep) { 7377f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 7387f1410a3SPeter Brune } 7397f1410a3SPeter Brune } 740bf7f4e0aSPeter Brune PetscFunctionReturn(0); 741bf7f4e0aSPeter Brune } 742bf7f4e0aSPeter Brune 743bf7f4e0aSPeter Brune #undef __FUNCT__ 744f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType" 745ea5d4fccSPeter Brune /*@C 746f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 747f40b411bSPeter Brune 748f40b411bSPeter Brune Input Parameters: 749f40b411bSPeter Brune + linesearch - linesearch context. 750f40b411bSPeter Brune - type - The type of line search to be used 751f40b411bSPeter Brune 752cd7522eaSPeter Brune Available Types: 753cd7522eaSPeter Brune + basic - Simple damping line search. 754cd7522eaSPeter Brune . bt - Backtracking line search over the L2 norm of the function. 755cd7522eaSPeter Brune . l2 - Secant line search over the L2 norm of the function. 756cd7522eaSPeter Brune . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x). 757cd7522eaSPeter Brune - shell - User provided SNESLineSearch implementation. 758cd7522eaSPeter Brune 759f1c6b773SPeter Brune Logically Collective on SNESLineSearch 760f40b411bSPeter Brune 761f40b411bSPeter Brune Level: intermediate 762f40b411bSPeter Brune 763f40b411bSPeter Brune 764f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 765f40b411bSPeter Brune @*/ 766f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type) 767bf7f4e0aSPeter Brune { 768bf7f4e0aSPeter Brune 769f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 770bf7f4e0aSPeter Brune PetscBool match; 771bf7f4e0aSPeter Brune 772bf7f4e0aSPeter Brune PetscFunctionBegin; 773f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 774bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 775bf7f4e0aSPeter Brune 776bf7f4e0aSPeter Brune ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 777bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 778bf7f4e0aSPeter Brune 779f1c6b773SPeter Brune ierr = PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 780bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 781bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 782bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 783bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 784bf7f4e0aSPeter Brune linesearch->ops->destroy = PETSC_NULL; 785bf7f4e0aSPeter Brune } 786f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 787bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 788bf7f4e0aSPeter Brune linesearch->ops->view = 0; 789bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 790bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 791bf7f4e0aSPeter Brune 792bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 793bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 794bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS) 795bf7f4e0aSPeter Brune if (PetscAMSPublishAll) { 796bf7f4e0aSPeter Brune ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 797bf7f4e0aSPeter Brune } 798bf7f4e0aSPeter Brune #endif 799bf7f4e0aSPeter Brune PetscFunctionReturn(0); 800bf7f4e0aSPeter Brune } 801bf7f4e0aSPeter Brune 802bf7f4e0aSPeter Brune #undef __FUNCT__ 803f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES" 804f40b411bSPeter Brune /*@ 805*78bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 806f40b411bSPeter Brune 807f40b411bSPeter Brune Input Parameters: 808f40b411bSPeter Brune + linesearch - linesearch context. 809f40b411bSPeter Brune - snes - The snes instance 810f40b411bSPeter Brune 811*78bcb3b5SPeter Brune Level: developer 812*78bcb3b5SPeter Brune 813*78bcb3b5SPeter Brune Notes: 814*78bcb3b5SPeter Brune This happens automatically when the line search is gotten/created with 815*78bcb3b5SPeter Brune SNESGetSNESLineSearch(). This routine is therefore mainly called within SNES 816*78bcb3b5SPeter Brune implementations. 817f40b411bSPeter Brune 818f40b411bSPeter Brune 819cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 820f40b411bSPeter Brune @*/ 821f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){ 822bf7f4e0aSPeter Brune PetscFunctionBegin; 823f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 824bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 825bf7f4e0aSPeter Brune linesearch->snes = snes; 826bf7f4e0aSPeter Brune PetscFunctionReturn(0); 827bf7f4e0aSPeter Brune } 828bf7f4e0aSPeter Brune 829bf7f4e0aSPeter Brune #undef __FUNCT__ 830f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES" 831f40b411bSPeter Brune /*@ 832*78bcb3b5SPeter Brune SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation. 833f40b411bSPeter Brune 834f40b411bSPeter Brune Input Parameters: 835f40b411bSPeter Brune . linesearch - linesearch context. 836f40b411bSPeter Brune 837f40b411bSPeter Brune Output Parameters: 838f40b411bSPeter Brune . snes - The snes instance 839f40b411bSPeter Brune 840*78bcb3b5SPeter Brune Level: advanced 841f40b411bSPeter Brune 842cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 843f40b411bSPeter Brune @*/ 844f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){ 845bf7f4e0aSPeter Brune PetscFunctionBegin; 846f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8476a388c36SPeter Brune PetscValidPointer(snes, 2); 848bf7f4e0aSPeter Brune *snes = linesearch->snes; 849bf7f4e0aSPeter Brune PetscFunctionReturn(0); 850bf7f4e0aSPeter Brune } 851bf7f4e0aSPeter Brune 8526a388c36SPeter Brune #undef __FUNCT__ 853f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda" 854f40b411bSPeter Brune /*@ 855f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 856f40b411bSPeter Brune 857f40b411bSPeter Brune Input Parameters: 858f40b411bSPeter Brune . linesearch - linesearch context. 859f40b411bSPeter Brune 860f40b411bSPeter Brune Output Parameters: 861cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 862f40b411bSPeter Brune 863*78bcb3b5SPeter Brune Level: advanced 864*78bcb3b5SPeter Brune 865*78bcb3b5SPeter Brune Notes: This is useful in methods where the solver is ill-scaled and 866*78bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 867*78bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 868*78bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 869*78bcb3b5SPeter Brune 870f40b411bSPeter Brune 871cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 872f40b411bSPeter Brune @*/ 873f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 8746a388c36SPeter Brune { 8756a388c36SPeter Brune PetscFunctionBegin; 876f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8776a388c36SPeter Brune PetscValidPointer(lambda, 2); 8786a388c36SPeter Brune *lambda = linesearch->lambda; 8796a388c36SPeter Brune PetscFunctionReturn(0); 8806a388c36SPeter Brune } 8816a388c36SPeter Brune 8826a388c36SPeter Brune #undef __FUNCT__ 883f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda" 884f40b411bSPeter Brune /*@ 885f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 886f40b411bSPeter Brune 887f40b411bSPeter Brune Input Parameters: 888f40b411bSPeter Brune + linesearch - linesearch context. 889f40b411bSPeter Brune - lambda - The last steplength. 890f40b411bSPeter Brune 891cd7522eaSPeter Brune Notes: 892cd7522eaSPeter Brune This routine is typically used within implementations of SNESLineSearchApply 893cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 894cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 895cd7522eaSPeter Brune as an inner scaling parameter. 896cd7522eaSPeter Brune 897*78bcb3b5SPeter Brune Level: advanced 898f40b411bSPeter Brune 899f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 900f40b411bSPeter Brune @*/ 901f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 9026a388c36SPeter Brune { 9036a388c36SPeter Brune PetscFunctionBegin; 904f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9056a388c36SPeter Brune linesearch->lambda = lambda; 9066a388c36SPeter Brune PetscFunctionReturn(0); 9076a388c36SPeter Brune } 9086a388c36SPeter Brune 9096a388c36SPeter Brune #undef __FUNCT__ 910f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances" 911f40b411bSPeter Brune /*@ 912*78bcb3b5SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the method. These include 913*78bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 914*78bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 915*78bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 916f40b411bSPeter Brune 917f40b411bSPeter Brune Input Parameters: 918f40b411bSPeter Brune . linesearch - linesearch context. 919f40b411bSPeter Brune 920f40b411bSPeter Brune Output Parameters: 921516fe3c3SPeter Brune + steptol - The minimum steplength 9226cc8e53bSPeter Brune . maxstep - The maximum steplength 923516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 924516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 925516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 926516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 927f40b411bSPeter Brune 928*78bcb3b5SPeter Brune Level: intermediate 929*78bcb3b5SPeter Brune 930*78bcb3b5SPeter Brune Notes: 931*78bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 932*78bcb3b5SPeter Brune the method requires. 933516fe3c3SPeter Brune 934f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 935f40b411bSPeter Brune @*/ 936f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 9376a388c36SPeter Brune { 9386a388c36SPeter Brune PetscFunctionBegin; 939f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 940516fe3c3SPeter Brune if (steptol) { 9416a388c36SPeter Brune PetscValidPointer(steptol, 2); 9426a388c36SPeter Brune *steptol = linesearch->steptol; 943516fe3c3SPeter Brune } 944516fe3c3SPeter Brune if (maxstep) { 945516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 946516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 947516fe3c3SPeter Brune } 948516fe3c3SPeter Brune if (rtol) { 949516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 950516fe3c3SPeter Brune *rtol = linesearch->rtol; 951516fe3c3SPeter Brune } 952516fe3c3SPeter Brune if (atol) { 953516fe3c3SPeter Brune PetscValidPointer(atol, 5); 954516fe3c3SPeter Brune *atol = linesearch->atol; 955516fe3c3SPeter Brune } 956516fe3c3SPeter Brune if (ltol) { 957516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 958516fe3c3SPeter Brune *ltol = linesearch->ltol; 959516fe3c3SPeter Brune } 960516fe3c3SPeter Brune if (max_its) { 961516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 962516fe3c3SPeter Brune *max_its = linesearch->max_its; 963516fe3c3SPeter Brune } 9646a388c36SPeter Brune PetscFunctionReturn(0); 9656a388c36SPeter Brune } 9666a388c36SPeter Brune 9676a388c36SPeter Brune #undef __FUNCT__ 968f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances" 969f40b411bSPeter Brune /*@ 970*78bcb3b5SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the method. These include 971*78bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 972*78bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 973*78bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 974f40b411bSPeter Brune 975f40b411bSPeter Brune Input Parameters: 976516fe3c3SPeter Brune + linesearch - linesearch context. 977516fe3c3SPeter Brune . steptol - The minimum steplength 9786cc8e53bSPeter Brune . maxstep - The maximum steplength 979516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 980516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 981516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 982516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 983f40b411bSPeter Brune 984*78bcb3b5SPeter Brune Notes: 985*78bcb3b5SPeter Brune The user may choose to not set any of the tolerances in the method using PETSC_DEFAULT in 986*78bcb3b5SPeter Brune place of an argument. 987f40b411bSPeter Brune 988*78bcb3b5SPeter Brune Level: intermediate 989516fe3c3SPeter Brune 990f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 991f40b411bSPeter Brune @*/ 992f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 9936a388c36SPeter Brune { 9946a388c36SPeter Brune PetscFunctionBegin; 995f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 996d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 997d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 998d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 999d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1000d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1001d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1002d3952184SSatish Balay 1003d3952184SSatish Balay if ( steptol!= PETSC_DEFAULT) { 1004d3952184SSatish Balay if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol); 10056a388c36SPeter Brune linesearch->steptol = steptol; 1006d3952184SSatish Balay } 1007d3952184SSatish Balay 1008d3952184SSatish Balay if ( maxstep!= PETSC_DEFAULT) { 1009d3952184SSatish Balay if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep); 1010516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1011d3952184SSatish Balay } 1012d3952184SSatish Balay 1013d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1014d3952184SSatish Balay if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol); 1015516fe3c3SPeter Brune linesearch->rtol = rtol; 1016d3952184SSatish Balay } 1017d3952184SSatish Balay 1018d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1019d3952184SSatish Balay if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol); 1020516fe3c3SPeter Brune linesearch->atol = atol; 1021d3952184SSatish Balay } 1022d3952184SSatish Balay 1023d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1024d3952184SSatish Balay if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol); 1025516fe3c3SPeter Brune linesearch->ltol = ltol; 1026d3952184SSatish Balay } 1027d3952184SSatish Balay 1028d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1029d3952184SSatish Balay if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1030516fe3c3SPeter Brune linesearch->max_its = max_its; 1031d3952184SSatish Balay } 1032d3952184SSatish Balay 10336a388c36SPeter Brune PetscFunctionReturn(0); 10346a388c36SPeter Brune } 10356a388c36SPeter Brune 1036516fe3c3SPeter Brune 10376a388c36SPeter Brune #undef __FUNCT__ 1038f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping" 1039f40b411bSPeter Brune /*@ 1040f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1041f40b411bSPeter Brune 1042f40b411bSPeter Brune Input Parameters: 1043f40b411bSPeter Brune . linesearch - linesearch context. 1044f40b411bSPeter Brune 1045f40b411bSPeter Brune Output Parameters: 1046f40b411bSPeter Brune . damping - The damping parameter. 1047f40b411bSPeter Brune 1048*78bcb3b5SPeter Brune Level: advanced 1049f40b411bSPeter Brune 1050*78bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1051f40b411bSPeter Brune @*/ 1052f40b411bSPeter Brune 1053f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 10546a388c36SPeter Brune { 10556a388c36SPeter Brune PetscFunctionBegin; 1056f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 10576a388c36SPeter Brune PetscValidPointer(damping, 2); 10586a388c36SPeter Brune *damping = linesearch->damping; 10596a388c36SPeter Brune PetscFunctionReturn(0); 10606a388c36SPeter Brune } 10616a388c36SPeter Brune 10626a388c36SPeter Brune #undef __FUNCT__ 1063f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping" 1064f40b411bSPeter Brune /*@ 1065f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1066f40b411bSPeter Brune 1067f40b411bSPeter Brune Input Parameters: 1068*78bcb3b5SPeter Brune . linesearch - linesearch context 1069*78bcb3b5SPeter Brune . damping - The damping parameter 1070f40b411bSPeter Brune 1071f40b411bSPeter Brune Level: intermediate 1072f40b411bSPeter Brune 1073cd7522eaSPeter Brune Notes: 1074cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1075cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 1076*78bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1077cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1078cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1079cd7522eaSPeter Brune 1080f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1081f40b411bSPeter Brune @*/ 1082f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 10836a388c36SPeter Brune { 10846a388c36SPeter Brune PetscFunctionBegin; 1085f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 10866a388c36SPeter Brune linesearch->damping = damping; 10876a388c36SPeter Brune PetscFunctionReturn(0); 10886a388c36SPeter Brune } 10896a388c36SPeter Brune 10906a388c36SPeter Brune #undef __FUNCT__ 109159405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder" 109259405d5eSPeter Brune /*@ 109359405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 109459405d5eSPeter Brune 109559405d5eSPeter Brune Input Parameters: 1096*78bcb3b5SPeter Brune . linesearch - linesearch context 109759405d5eSPeter Brune 109859405d5eSPeter Brune Output Parameters: 109959405d5eSPeter Brune . order - The order. 110059405d5eSPeter Brune 1101*78bcb3b5SPeter Brune Possible Values for order: 1102*78bcb3b5SPeter Brune + SNES_LINESEARCH_LINEAR 1103*78bcb3b5SPeter Brune . SNES_LINESEARCH_CUBIC 1104*78bcb3b5SPeter Brune - SNES_LINESEARCH_QUADRATIC 1105*78bcb3b5SPeter Brune 110659405d5eSPeter Brune Level: intermediate 110759405d5eSPeter Brune 110859405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 110959405d5eSPeter Brune @*/ 111059405d5eSPeter Brune 111159405d5eSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,SNESLineSearchOrder *order) 111259405d5eSPeter Brune { 111359405d5eSPeter Brune PetscFunctionBegin; 111459405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 111559405d5eSPeter Brune PetscValidPointer(order, 2); 111659405d5eSPeter Brune *order = linesearch->order; 111759405d5eSPeter Brune PetscFunctionReturn(0); 111859405d5eSPeter Brune } 111959405d5eSPeter Brune 112059405d5eSPeter Brune #undef __FUNCT__ 112159405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder" 112259405d5eSPeter Brune /*@ 112359405d5eSPeter Brune SNESLineSearchSetOrder - Sets the line search damping paramter. 112459405d5eSPeter Brune 112559405d5eSPeter Brune Input Parameters: 1126*78bcb3b5SPeter Brune . linesearch - linesearch context 1127*78bcb3b5SPeter Brune . order - The damping parameter 112859405d5eSPeter Brune 112959405d5eSPeter Brune Level: intermediate 113059405d5eSPeter Brune 1131*78bcb3b5SPeter Brune Possible Values for order: 1132*78bcb3b5SPeter Brune + SNES_LINESEARCH_LINEAR 1133*78bcb3b5SPeter Brune . SNES_LINESEARCH_CUBIC 1134*78bcb3b5SPeter Brune - SNES_LINESEARCH_QUADRATIC 1135*78bcb3b5SPeter Brune 113659405d5eSPeter Brune Notes: 113759405d5eSPeter Brune Variable orders are supported by the following line searches: 1138*78bcb3b5SPeter Brune + bt - cubic and quadratic 1139*78bcb3b5SPeter Brune - cp - linear and quadratic 114059405d5eSPeter Brune 114159405d5eSPeter Brune .seealso: SNESLineSearchGetOrder() 114259405d5eSPeter Brune @*/ 114359405d5eSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,SNESLineSearchOrder order) 114459405d5eSPeter Brune { 114559405d5eSPeter Brune PetscFunctionBegin; 114659405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 114759405d5eSPeter Brune linesearch->order = order; 114859405d5eSPeter Brune PetscFunctionReturn(0); 114959405d5eSPeter Brune } 115059405d5eSPeter Brune 115159405d5eSPeter Brune #undef __FUNCT__ 1152f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms" 1153f40b411bSPeter Brune /*@ 1154f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1155f40b411bSPeter Brune 1156f40b411bSPeter Brune Input Parameters: 1157*78bcb3b5SPeter Brune . linesearch - linesearch context 1158f40b411bSPeter Brune 1159f40b411bSPeter Brune Output Parameters: 1160f40b411bSPeter Brune + xnorm - The norm of the current solution 1161f40b411bSPeter Brune . fnorm - The norm of the current function 1162f40b411bSPeter Brune - ynorm - The norm of the current update 1163f40b411bSPeter Brune 1164cd7522eaSPeter Brune Notes: 1165cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1166cd7522eaSPeter Brune 1167*78bcb3b5SPeter Brune Level: developer 1168f40b411bSPeter Brune 1169f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1170f40b411bSPeter Brune @*/ 1171f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1172bf7f4e0aSPeter Brune { 1173bf7f4e0aSPeter Brune PetscFunctionBegin; 1174f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1175bf7f4e0aSPeter Brune if (xnorm) { 1176bf7f4e0aSPeter Brune *xnorm = linesearch->xnorm; 1177bf7f4e0aSPeter Brune } 1178bf7f4e0aSPeter Brune if (fnorm) { 1179bf7f4e0aSPeter Brune *fnorm = linesearch->fnorm; 1180bf7f4e0aSPeter Brune } 1181bf7f4e0aSPeter Brune if (ynorm) { 1182bf7f4e0aSPeter Brune *ynorm = linesearch->ynorm; 1183bf7f4e0aSPeter Brune } 1184bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1185bf7f4e0aSPeter Brune } 1186bf7f4e0aSPeter Brune 1187e7058c64SPeter Brune #undef __FUNCT__ 1188f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms" 1189f40b411bSPeter Brune /*@ 1190f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1191f40b411bSPeter Brune 1192f40b411bSPeter Brune Input Parameters: 1193*78bcb3b5SPeter Brune + linesearch - linesearch context 1194f40b411bSPeter Brune . xnorm - The norm of the current solution 1195f40b411bSPeter Brune . fnorm - The norm of the current function 1196f40b411bSPeter Brune - ynorm - The norm of the current update 1197f40b411bSPeter Brune 1198*78bcb3b5SPeter Brune Level: advanced 1199f40b411bSPeter Brune 1200f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1201f40b411bSPeter Brune @*/ 1202f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 12036a388c36SPeter Brune { 12046a388c36SPeter Brune PetscFunctionBegin; 1205f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12066a388c36SPeter Brune linesearch->xnorm = xnorm; 12076a388c36SPeter Brune linesearch->fnorm = fnorm; 12086a388c36SPeter Brune linesearch->ynorm = ynorm; 12096a388c36SPeter Brune PetscFunctionReturn(0); 12106a388c36SPeter Brune } 12116a388c36SPeter Brune 12126a388c36SPeter Brune #undef __FUNCT__ 1213f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms" 1214f40b411bSPeter Brune /*@ 1215f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1216f40b411bSPeter Brune 1217f40b411bSPeter Brune Input Parameters: 1218*78bcb3b5SPeter Brune . linesearch - linesearch context 1219f40b411bSPeter Brune 1220f40b411bSPeter Brune Options Database Keys: 1221cd7522eaSPeter Brune . -snes_linesearch_norms - turn norm computation on or off. 1222f40b411bSPeter Brune 1223f40b411bSPeter Brune Level: intermediate 1224f40b411bSPeter Brune 1225*78bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1226f40b411bSPeter Brune @*/ 1227f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 12286a388c36SPeter Brune { 12296a388c36SPeter Brune PetscErrorCode ierr; 12309bd66eb0SPeter Brune SNES snes; 12316a388c36SPeter Brune PetscFunctionBegin; 12326a388c36SPeter Brune if (linesearch->norms) { 12339bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1234f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 12359bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12369bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12379bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 12389bd66eb0SPeter Brune } else { 12396a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12406a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12416a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12426a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12436a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12446a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12456a388c36SPeter Brune } 12469bd66eb0SPeter Brune } 12476a388c36SPeter Brune PetscFunctionReturn(0); 12486a388c36SPeter Brune } 12496a388c36SPeter Brune 12506f263ca3SPeter Brune 12516f263ca3SPeter Brune #undef __FUNCT__ 12526f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms" 12536f263ca3SPeter Brune /*@ 12546f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 12556f263ca3SPeter Brune 12566f263ca3SPeter Brune Input Parameters: 1257*78bcb3b5SPeter Brune + linesearch - linesearch context 1258*78bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 12596f263ca3SPeter Brune 12606f263ca3SPeter Brune Options Database Keys: 12616f263ca3SPeter Brune . -snes_linesearch_norms - turn norm computation on or off. 12626f263ca3SPeter Brune 12636f263ca3SPeter Brune Notes: 12646f263ca3SPeter Brune This is most relevant to the SNES_LINESEARCH_BASIC line search type. 12656f263ca3SPeter Brune 12666f263ca3SPeter Brune Level: intermediate 12676f263ca3SPeter Brune 12686f263ca3SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNES_LINESEARCH_BASIC 12696f263ca3SPeter Brune @*/ 12706f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 12716f263ca3SPeter Brune { 12726f263ca3SPeter Brune PetscFunctionBegin; 12736f263ca3SPeter Brune linesearch->norms = flg; 12746f263ca3SPeter Brune PetscFunctionReturn(0); 12756f263ca3SPeter Brune } 12766f263ca3SPeter Brune 12776a388c36SPeter Brune #undef __FUNCT__ 1278f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs" 1279f40b411bSPeter Brune /*@ 1280f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1281f40b411bSPeter Brune 1282f40b411bSPeter Brune Input Parameters: 1283*78bcb3b5SPeter Brune . linesearch - linesearch context 1284f40b411bSPeter Brune 1285f40b411bSPeter Brune Output Parameters: 1286f40b411bSPeter Brune + X - The old solution 1287f40b411bSPeter Brune . F - The old function 1288f40b411bSPeter Brune . Y - The search direction 1289f40b411bSPeter Brune . W - The new solution 1290f40b411bSPeter Brune - G - The new function 1291f40b411bSPeter Brune 1292*78bcb3b5SPeter Brune Level: advanced 1293f40b411bSPeter Brune 1294f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1295f40b411bSPeter Brune @*/ 1296f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) { 12976a388c36SPeter Brune PetscFunctionBegin; 1298f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12996a388c36SPeter Brune if (X) { 13006a388c36SPeter Brune PetscValidPointer(X, 2); 13016a388c36SPeter Brune *X = linesearch->vec_sol; 13026a388c36SPeter Brune } 13036a388c36SPeter Brune if (F) { 13046a388c36SPeter Brune PetscValidPointer(F, 3); 13056a388c36SPeter Brune *F = linesearch->vec_func; 13066a388c36SPeter Brune } 13076a388c36SPeter Brune if (Y) { 13086a388c36SPeter Brune PetscValidPointer(Y, 4); 13096a388c36SPeter Brune *Y = linesearch->vec_update; 13106a388c36SPeter Brune } 13116a388c36SPeter Brune if (W) { 13126a388c36SPeter Brune PetscValidPointer(W, 5); 13136a388c36SPeter Brune *W = linesearch->vec_sol_new; 13146a388c36SPeter Brune } 13156a388c36SPeter Brune if (G) { 13166a388c36SPeter Brune PetscValidPointer(G, 6); 13176a388c36SPeter Brune *G = linesearch->vec_func_new; 13186a388c36SPeter Brune } 13196a388c36SPeter Brune 13206a388c36SPeter Brune PetscFunctionReturn(0); 13216a388c36SPeter Brune } 13226a388c36SPeter Brune 13236a388c36SPeter Brune #undef __FUNCT__ 1324f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs" 1325f40b411bSPeter Brune /*@ 1326f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1327f40b411bSPeter Brune 1328f40b411bSPeter Brune Input Parameters: 1329*78bcb3b5SPeter Brune + linesearch - linesearch context 1330f40b411bSPeter Brune . X - The old solution 1331f40b411bSPeter Brune . F - The old function 1332f40b411bSPeter Brune . Y - The search direction 1333f40b411bSPeter Brune . W - The new solution 1334f40b411bSPeter Brune - G - The new function 1335f40b411bSPeter Brune 1336*78bcb3b5SPeter Brune Level: advanced 1337f40b411bSPeter Brune 1338f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1339f40b411bSPeter Brune @*/ 1340f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) { 13416a388c36SPeter Brune PetscFunctionBegin; 1342f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13436a388c36SPeter Brune if (X) { 13446a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 13456a388c36SPeter Brune linesearch->vec_sol = X; 13466a388c36SPeter Brune } 13476a388c36SPeter Brune if (F) { 13486a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 13496a388c36SPeter Brune linesearch->vec_func = F; 13506a388c36SPeter Brune } 13516a388c36SPeter Brune if (Y) { 13526a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 13536a388c36SPeter Brune linesearch->vec_update = Y; 13546a388c36SPeter Brune } 13556a388c36SPeter Brune if (W) { 13566a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 13576a388c36SPeter Brune linesearch->vec_sol_new = W; 13586a388c36SPeter Brune } 13596a388c36SPeter Brune if (G) { 13606a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 13616a388c36SPeter Brune linesearch->vec_func_new = G; 13626a388c36SPeter Brune } 13636a388c36SPeter Brune 13646a388c36SPeter Brune PetscFunctionReturn(0); 13656a388c36SPeter Brune } 13666a388c36SPeter Brune 13676a388c36SPeter Brune #undef __FUNCT__ 1368f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1369e7058c64SPeter Brune /*@C 1370f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1371e7058c64SPeter Brune SNES options in the database. 1372e7058c64SPeter Brune 1373cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1374e7058c64SPeter Brune 1375e7058c64SPeter Brune Input Parameters: 1376e7058c64SPeter Brune + snes - the SNES context 1377e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1378e7058c64SPeter Brune 1379e7058c64SPeter Brune Notes: 1380e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1381e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1382e7058c64SPeter Brune 1383e7058c64SPeter Brune Level: advanced 1384e7058c64SPeter Brune 1385f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1386e7058c64SPeter Brune 1387e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1388e7058c64SPeter Brune @*/ 1389f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1390e7058c64SPeter Brune { 1391e7058c64SPeter Brune PetscErrorCode ierr; 1392e7058c64SPeter Brune 1393e7058c64SPeter Brune PetscFunctionBegin; 1394f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1395e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1396e7058c64SPeter Brune PetscFunctionReturn(0); 1397e7058c64SPeter Brune } 1398e7058c64SPeter Brune 1399e7058c64SPeter Brune #undef __FUNCT__ 1400f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1401e7058c64SPeter Brune /*@C 1402f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1403f1c6b773SPeter Brune SNESLineSearch options in the database. 1404e7058c64SPeter Brune 1405e7058c64SPeter Brune Not Collective 1406e7058c64SPeter Brune 1407e7058c64SPeter Brune Input Parameter: 1408cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1409e7058c64SPeter Brune 1410e7058c64SPeter Brune Output Parameter: 1411e7058c64SPeter Brune . prefix - pointer to the prefix string used 1412e7058c64SPeter Brune 1413e7058c64SPeter Brune Notes: On the fortran side, the user should pass in a string 'prefix' of 1414e7058c64SPeter Brune sufficient length to hold the prefix. 1415e7058c64SPeter Brune 1416e7058c64SPeter Brune Level: advanced 1417e7058c64SPeter Brune 1418f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1419e7058c64SPeter Brune 1420e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1421e7058c64SPeter Brune @*/ 1422f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1423e7058c64SPeter Brune { 1424e7058c64SPeter Brune PetscErrorCode ierr; 1425e7058c64SPeter Brune 1426e7058c64SPeter Brune PetscFunctionBegin; 1427f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1428e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1429e7058c64SPeter Brune PetscFunctionReturn(0); 1430e7058c64SPeter Brune } 1431bf7f4e0aSPeter Brune 1432bf7f4e0aSPeter Brune #undef __FUNCT__ 1433f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork" 1434f40b411bSPeter Brune /*@ 1435f1c6b773SPeter Brune SNESLineSearchGetWork - Gets work vectors for the line search. 1436f40b411bSPeter Brune 1437f40b411bSPeter Brune Input Parameter: 1438f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1439f40b411bSPeter Brune - nwork - the number of work vectors 1440f40b411bSPeter Brune 1441f40b411bSPeter Brune Level: developer 1442f40b411bSPeter Brune 1443cd7522eaSPeter Brune Notes: 1444cd7522eaSPeter Brune This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation. 1445cd7522eaSPeter Brune 1446f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1447f40b411bSPeter Brune 1448f40b411bSPeter Brune .seealso: SNESDefaultGetWork() 1449f40b411bSPeter Brune @*/ 1450f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork) 1451bf7f4e0aSPeter Brune { 1452bf7f4e0aSPeter Brune PetscErrorCode ierr; 1453bf7f4e0aSPeter Brune PetscFunctionBegin; 1454bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1455bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1456bf7f4e0aSPeter Brune } else { 1457bf7f4e0aSPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1458bf7f4e0aSPeter Brune } 1459bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1460bf7f4e0aSPeter Brune } 1461bf7f4e0aSPeter Brune 14626a388c36SPeter Brune 1463bf7f4e0aSPeter Brune #undef __FUNCT__ 1464f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess" 1465f40b411bSPeter Brune /*@ 1466f1c6b773SPeter Brune SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application 1467f40b411bSPeter Brune 1468f40b411bSPeter Brune Input Parameters: 1469*78bcb3b5SPeter Brune . linesearch - linesearch context 1470f40b411bSPeter Brune 1471f40b411bSPeter Brune Output Parameters: 1472f40b411bSPeter Brune . success - The success or failure status. 1473f40b411bSPeter Brune 1474cd7522eaSPeter Brune Notes: 1475cd7522eaSPeter Brune This is typically called after SNESLineSearchApply in order to determine if the line-search failed 1476cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1477cd7522eaSPeter Brune 1478f40b411bSPeter Brune Level: intermediate 1479f40b411bSPeter Brune 1480f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess() 1481f40b411bSPeter Brune @*/ 1482f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success) 1483bf7f4e0aSPeter Brune { 1484bf7f4e0aSPeter Brune PetscFunctionBegin; 1485f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14866a388c36SPeter Brune PetscValidPointer(success, 2); 1487bf7f4e0aSPeter Brune if (success) { 1488bf7f4e0aSPeter Brune *success = linesearch->success; 1489bf7f4e0aSPeter Brune } 1490bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1491bf7f4e0aSPeter Brune } 1492bf7f4e0aSPeter Brune 1493bf7f4e0aSPeter Brune #undef __FUNCT__ 1494f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess" 1495f40b411bSPeter Brune /*@ 1496f1c6b773SPeter Brune SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application 1497f40b411bSPeter Brune 1498f40b411bSPeter Brune Input Parameters: 1499*78bcb3b5SPeter Brune + linesearch - linesearch context 1500f40b411bSPeter Brune - success - The success or failure status. 1501f40b411bSPeter Brune 1502cd7522eaSPeter Brune Notes: 1503cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1504cd7522eaSPeter Brune the success or failure of the line search method. 1505cd7522eaSPeter Brune 1506*78bcb3b5SPeter Brune Level: developer 1507f40b411bSPeter Brune 1508f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess() 1509f40b411bSPeter Brune @*/ 1510f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success) 15116a388c36SPeter Brune { 1512f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15136a388c36SPeter Brune PetscFunctionBegin; 15146a388c36SPeter Brune linesearch->success = success; 15156a388c36SPeter Brune PetscFunctionReturn(0); 15166a388c36SPeter Brune } 15176a388c36SPeter Brune 15186a388c36SPeter Brune #undef __FUNCT__ 1519f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions" 15209bd66eb0SPeter Brune /*@C 1521f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 15229bd66eb0SPeter Brune 15239bd66eb0SPeter Brune Input Parameters: 15249bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 15259bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 15269bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 15279bd66eb0SPeter Brune 15289bd66eb0SPeter Brune Logically Collective on SNES 15299bd66eb0SPeter Brune 15309bd66eb0SPeter Brune Calling sequence of projectfunc: 15319bd66eb0SPeter Brune .vb 15329bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 15339bd66eb0SPeter Brune .ve 15349bd66eb0SPeter Brune 15359bd66eb0SPeter Brune Input parameters for projectfunc: 15369bd66eb0SPeter Brune + snes - nonlinear context 15379bd66eb0SPeter Brune - X - current solution 15389bd66eb0SPeter Brune 1539cd7522eaSPeter Brune Output parameters for projectfunc: 15409bd66eb0SPeter Brune . X - Projected solution 15419bd66eb0SPeter Brune 15429bd66eb0SPeter Brune Calling sequence of normfunc: 15439bd66eb0SPeter Brune .vb 15449bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 15459bd66eb0SPeter Brune .ve 15469bd66eb0SPeter Brune 1547cd7522eaSPeter Brune Input parameters for normfunc: 15489bd66eb0SPeter Brune + snes - nonlinear context 15499bd66eb0SPeter Brune . X - current solution 15509bd66eb0SPeter Brune - F - current residual 15519bd66eb0SPeter Brune 1552cd7522eaSPeter Brune Output parameters for normfunc: 15539bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 15549bd66eb0SPeter Brune 1555cd7522eaSPeter Brune Notes: 1556cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1557cd7522eaSPeter Brune 1558cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1559cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 15609bd66eb0SPeter Brune 15619bd66eb0SPeter Brune Level: developer 15629bd66eb0SPeter Brune 15639bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 15649bd66eb0SPeter Brune 1565f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 15669bd66eb0SPeter Brune @*/ 1567f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 15689bd66eb0SPeter Brune { 15699bd66eb0SPeter Brune PetscFunctionBegin; 1570f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15719bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 15729bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 15739bd66eb0SPeter Brune PetscFunctionReturn(0); 15749bd66eb0SPeter Brune } 15759bd66eb0SPeter Brune 15769bd66eb0SPeter Brune #undef __FUNCT__ 1577f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions" 15789bd66eb0SPeter Brune /*@C 1579f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 15809bd66eb0SPeter Brune 15819bd66eb0SPeter Brune Input Parameters: 15829bd66eb0SPeter Brune . snes - nonlinear context obtained from SNESCreate() 15839bd66eb0SPeter Brune 15849bd66eb0SPeter Brune Output Parameters: 15859bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 15869bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 15879bd66eb0SPeter Brune 15889bd66eb0SPeter Brune Logically Collective on SNES 15899bd66eb0SPeter Brune 15909bd66eb0SPeter Brune Level: developer 15919bd66eb0SPeter Brune 15929bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 15939bd66eb0SPeter Brune 1594f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 15959bd66eb0SPeter Brune @*/ 1596f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 15979bd66eb0SPeter Brune { 15989bd66eb0SPeter Brune PetscFunctionBegin; 15999bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16009bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16019bd66eb0SPeter Brune PetscFunctionReturn(0); 16029bd66eb0SPeter Brune } 16039bd66eb0SPeter Brune 16049bd66eb0SPeter Brune #undef __FUNCT__ 1605f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister" 1606bf7f4e0aSPeter Brune /*@C 1607f1c6b773SPeter Brune SNESLineSearchRegister - See SNESLineSearchRegisterDynamic() 1608bf7f4e0aSPeter Brune 1609bf7f4e0aSPeter Brune Level: advanced 1610bf7f4e0aSPeter Brune @*/ 1611f1c6b773SPeter Brune PetscErrorCode SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch)) 1612bf7f4e0aSPeter Brune { 1613bf7f4e0aSPeter Brune char fullname[PETSC_MAX_PATH_LEN]; 1614bf7f4e0aSPeter Brune PetscErrorCode ierr; 1615bf7f4e0aSPeter Brune 1616bf7f4e0aSPeter Brune PetscFunctionBegin; 1617bf7f4e0aSPeter Brune ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1618f1c6b773SPeter Brune ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1619bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1620bf7f4e0aSPeter Brune } 1621