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: 208e557f58SPeter Brune . outlinesearch - the new linesearch context 21f40b411bSPeter Brune 2278bcb3b5SPeter Brune Level: beginner 23f40b411bSPeter Brune 24cd7522eaSPeter Brune .keywords: LineSearch, create, context 25f40b411bSPeter Brune 26f40b411bSPeter Brune .seealso: LineSearchDestroy() 27f40b411bSPeter Brune @*/ 28f40b411bSPeter Brune 29*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 30*bf388a1fSBarry Smith { 31bf7f4e0aSPeter Brune PetscErrorCode ierr; 32f1c6b773SPeter Brune SNESLineSearch linesearch; 33*bf388a1fSBarry Smith 34bf7f4e0aSPeter Brune PetscFunctionBegin; 35ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 36ea5d4fccSPeter Brune *outlinesearch = PETSC_NULL; 37*bf388a1fSBarry Smith ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, 0,"SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 38bf7f4e0aSPeter Brune 39bf7f4e0aSPeter Brune linesearch->ops->precheckstep = PETSC_NULL; 40bf7f4e0aSPeter Brune linesearch->ops->postcheckstep = PETSC_NULL; 41bf7f4e0aSPeter Brune 429bd66eb0SPeter Brune linesearch->vec_sol_new = PETSC_NULL; 439bd66eb0SPeter Brune linesearch->vec_func_new = PETSC_NULL; 449bd66eb0SPeter Brune linesearch->vec_sol = PETSC_NULL; 459bd66eb0SPeter Brune linesearch->vec_func = PETSC_NULL; 469bd66eb0SPeter Brune linesearch->vec_update = PETSC_NULL; 479bd66eb0SPeter Brune 48bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 49bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 50bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 51bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 52bf7f4e0aSPeter Brune linesearch->success = PETSC_TRUE; 53bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 54bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 55bf7f4e0aSPeter Brune linesearch->damping = 1.0; 56bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 57bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 58516fe3c3SPeter Brune linesearch->rtol = 1e-8; 59516fe3c3SPeter Brune linesearch->atol = 1e-15; 60516fe3c3SPeter Brune linesearch->ltol = 1e-8; 61bf7f4e0aSPeter Brune linesearch->precheckctx = PETSC_NULL; 62bf7f4e0aSPeter Brune linesearch->postcheckctx = PETSC_NULL; 6359405d5eSPeter Brune linesearch->max_its = 1; 64bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 65bf7f4e0aSPeter Brune *outlinesearch = linesearch; 66bf7f4e0aSPeter Brune PetscFunctionReturn(0); 67bf7f4e0aSPeter Brune } 68bf7f4e0aSPeter Brune 69bf7f4e0aSPeter Brune #undef __FUNCT__ 70f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp" 71f40b411bSPeter Brune /*@ 7278bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 7378bcb3b5SPeter Brune any required vectors. 74f40b411bSPeter Brune 75cd7522eaSPeter Brune Collective on SNESLineSearch 76f40b411bSPeter Brune 77f40b411bSPeter Brune Input Parameters: 78f40b411bSPeter Brune . linesearch - The LineSearch instance. 79f40b411bSPeter Brune 80cd7522eaSPeter Brune Notes: 81cd7522eaSPeter Brune For most cases, this needn't be called outside of SNESLineSearchApply(). 82cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 8378bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 84cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 85cd7522eaSPeter Brune allocated upfront. 86cd7522eaSPeter Brune 8778bcb3b5SPeter Brune Level: advanced 88f40b411bSPeter Brune 89f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 90f40b411bSPeter Brune 91f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 92f40b411bSPeter Brune @*/ 93f40b411bSPeter Brune 94*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 95*bf388a1fSBarry Smith { 96bf7f4e0aSPeter Brune PetscErrorCode ierr; 97*bf388a1fSBarry Smith 98bf7f4e0aSPeter Brune PetscFunctionBegin; 99bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 1001a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 101bf7f4e0aSPeter Brune } 102bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 1039bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 104bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 1059bd66eb0SPeter Brune } 1069bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 1079bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 1089bd66eb0SPeter Brune } 109bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 110bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 111bf7f4e0aSPeter Brune } 112bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 113bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 114bf7f4e0aSPeter Brune } 115bf7f4e0aSPeter Brune PetscFunctionReturn(0); 116bf7f4e0aSPeter Brune } 117bf7f4e0aSPeter Brune 118bf7f4e0aSPeter Brune #undef __FUNCT__ 119f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset" 120f40b411bSPeter Brune 121f40b411bSPeter Brune /*@ 12278bcb3b5SPeter Brune SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search. 123f40b411bSPeter Brune 124f1c6b773SPeter Brune Collective on SNESLineSearch 125f40b411bSPeter Brune 126f40b411bSPeter Brune Input Parameters: 127f40b411bSPeter Brune . linesearch - The LineSearch instance. 128f40b411bSPeter Brune 12978bcb3b5SPeter Brune Level: intermediate 130f40b411bSPeter Brune 131cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 132f40b411bSPeter Brune 133f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 134f40b411bSPeter Brune @*/ 135f40b411bSPeter Brune 136*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 137*bf388a1fSBarry Smith { 138bf7f4e0aSPeter Brune PetscErrorCode ierr; 139*bf388a1fSBarry Smith 140bf7f4e0aSPeter Brune PetscFunctionBegin; 141bf7f4e0aSPeter Brune if (linesearch->ops->reset) { 142bf7f4e0aSPeter Brune (*linesearch->ops->reset)(linesearch); 143bf7f4e0aSPeter Brune } 144bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 145bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 146bf7f4e0aSPeter Brune 147bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 148bf7f4e0aSPeter Brune linesearch->nwork = 0; 149bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 150bf7f4e0aSPeter Brune PetscFunctionReturn(0); 151bf7f4e0aSPeter Brune } 152bf7f4e0aSPeter Brune 15386d74e61SPeter Brune #undef __FUNCT__ 154f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck" 15586d74e61SPeter Brune /*@C 156f1c6b773SPeter Brune SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine. 15786d74e61SPeter Brune 158f1c6b773SPeter Brune Logically Collective on SNESLineSearch 15986d74e61SPeter Brune 16086d74e61SPeter Brune Input Parameters: 161f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 16286d74e61SPeter Brune . func - [optional] function evaluation routine 16386d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 16486d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 16586d74e61SPeter Brune 16686d74e61SPeter Brune Calling sequence of func: 167f1c6b773SPeter Brune $ func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 16886d74e61SPeter Brune 16986d74e61SPeter Brune + x - solution vector 17086d74e61SPeter Brune . y - search direction vector 171cd7522eaSPeter Brune - changed - flag to indicate the precheck changed x or y. 17286d74e61SPeter Brune 17386d74e61SPeter Brune Level: intermediate 17486d74e61SPeter Brune 17586d74e61SPeter Brune .keywords: set, linesearch, pre-check 17686d74e61SPeter Brune 177f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck() 17886d74e61SPeter Brune @*/ 179f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx) 18086d74e61SPeter Brune { 1819bd66eb0SPeter Brune PetscFunctionBegin; 182f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 18386d74e61SPeter Brune if (func) linesearch->ops->precheckstep = func; 18486d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 18586d74e61SPeter Brune PetscFunctionReturn(0); 18686d74e61SPeter Brune } 18786d74e61SPeter Brune 18886d74e61SPeter Brune #undef __FUNCT__ 189f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck" 19086d74e61SPeter Brune /*@C 191cd7522eaSPeter Brune SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine. 19286d74e61SPeter Brune 19386d74e61SPeter Brune Input Parameters: 194f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 19586d74e61SPeter Brune 19686d74e61SPeter Brune Output Parameters: 19786d74e61SPeter Brune + func - [optional] function evaluation routine 19886d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 19986d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 20086d74e61SPeter Brune 20186d74e61SPeter Brune Level: intermediate 20286d74e61SPeter Brune 20386d74e61SPeter Brune .keywords: get, linesearch, pre-check 20486d74e61SPeter Brune 205f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 20686d74e61SPeter Brune @*/ 207f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx) 20886d74e61SPeter Brune { 2099bd66eb0SPeter Brune PetscFunctionBegin; 210f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 21186d74e61SPeter Brune if (func) *func = linesearch->ops->precheckstep; 21286d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 21386d74e61SPeter Brune PetscFunctionReturn(0); 21486d74e61SPeter Brune } 21586d74e61SPeter Brune 21686d74e61SPeter Brune #undef __FUNCT__ 217f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck" 21886d74e61SPeter Brune /*@C 219f1c6b773SPeter Brune SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine. 22086d74e61SPeter Brune 221f1c6b773SPeter Brune Logically Collective on SNESLineSearch 22286d74e61SPeter Brune 22386d74e61SPeter Brune Input Parameters: 224f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 22586d74e61SPeter Brune . func - [optional] function evaluation routine 22686d74e61SPeter Brune - ctx - [optional] user-defined context for private data for the 22786d74e61SPeter Brune function evaluation routine (may be PETSC_NULL) 22886d74e61SPeter Brune 22986d74e61SPeter Brune Calling sequence of func: 230fd573748SDmitry Karpeev $ func (SNESLineSearch linesearch,Vec x,Vec y, Vec w, *changed_y, PetscBool *changed_w); 23186d74e61SPeter Brune 23286d74e61SPeter Brune + x - old solution vector 23386d74e61SPeter Brune . y - search direction vector 23486d74e61SPeter Brune . w - new solution vector 2358e557f58SPeter Brune . changed_y - indicates that the line search changed y 2368e557f58SPeter Brune . changed_w - indicates that the line search changed w 23786d74e61SPeter Brune 23886d74e61SPeter Brune Level: intermediate 23986d74e61SPeter Brune 24086d74e61SPeter Brune .keywords: set, linesearch, post-check 24186d74e61SPeter Brune 242f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 24386d74e61SPeter Brune @*/ 244f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx) 24586d74e61SPeter Brune { 24686d74e61SPeter Brune PetscFunctionBegin; 247f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 24886d74e61SPeter Brune if (func) linesearch->ops->postcheckstep = func; 24986d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 25086d74e61SPeter Brune PetscFunctionReturn(0); 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 281bf7f4e0aSPeter Brune #undef __FUNCT__ 282f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck" 283f40b411bSPeter Brune /*@ 284f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 285f40b411bSPeter Brune 286cd7522eaSPeter Brune Logically Collective on SNESLineSearch 287f40b411bSPeter Brune 288f40b411bSPeter Brune Input Parameters: 2897b1df9c1SPeter Brune + linesearch - The linesearch instance. 2907b1df9c1SPeter Brune . X - The current solution 2917b1df9c1SPeter Brune - Y - The step direction 292f40b411bSPeter Brune 293f40b411bSPeter Brune Output Parameters: 2948e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 295f40b411bSPeter Brune 296f40b411bSPeter Brune Level: Beginner 297f40b411bSPeter Brune 298f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 299f40b411bSPeter Brune 300f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 301f40b411bSPeter Brune @*/ 3027b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 303bf7f4e0aSPeter Brune { 304bf7f4e0aSPeter Brune PetscErrorCode ierr; 305bf7f4e0aSPeter Brune PetscFunctionBegin; 306bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 307bf7f4e0aSPeter Brune if (linesearch->ops->precheckstep) { 3087b1df9c1SPeter Brune ierr = (*linesearch->ops->precheckstep)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr); 30938bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 310bf7f4e0aSPeter Brune } 311bf7f4e0aSPeter Brune PetscFunctionReturn(0); 312bf7f4e0aSPeter Brune } 313bf7f4e0aSPeter Brune 314bf7f4e0aSPeter Brune #undef __FUNCT__ 315f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck" 316f40b411bSPeter Brune /*@ 317f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 318f40b411bSPeter Brune 319cd7522eaSPeter Brune Logically Collective on SNESLineSearch 320f40b411bSPeter Brune 321f40b411bSPeter Brune Input Parameters: 3227b1df9c1SPeter Brune + linesearch - The linesearch context 3237b1df9c1SPeter Brune . X - The last solution 3247b1df9c1SPeter Brune . Y - The step direction 3257b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 326f40b411bSPeter Brune 327f40b411bSPeter Brune Output Parameters: 32878bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 32978bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 330f40b411bSPeter Brune 331f40b411bSPeter Brune Level: Intermediate 332f40b411bSPeter Brune 333f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 334f40b411bSPeter Brune 335f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 336f40b411bSPeter Brune @*/ 3377b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 338bf7f4e0aSPeter Brune { 339bf7f4e0aSPeter Brune PetscErrorCode ierr; 340*bf388a1fSBarry Smith 341bf7f4e0aSPeter Brune PetscFunctionBegin; 342bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 343bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 344bf7f4e0aSPeter Brune if (linesearch->ops->postcheckstep) { 3457b1df9c1SPeter Brune ierr = (*linesearch->ops->postcheckstep)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr); 34638bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 34738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 34886d74e61SPeter Brune } 34986d74e61SPeter Brune PetscFunctionReturn(0); 35086d74e61SPeter Brune } 35186d74e61SPeter Brune 35286d74e61SPeter Brune #undef __FUNCT__ 353f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard" 35486d74e61SPeter Brune /*@C 35586d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 35686d74e61SPeter Brune 357cd7522eaSPeter Brune Logically Collective on SNESLineSearch 35886d74e61SPeter Brune 35986d74e61SPeter Brune Input Arguments: 36086d74e61SPeter Brune + linesearch - linesearch context 36186d74e61SPeter Brune . X - base state for this step 36286d74e61SPeter Brune . Y - initial correction 36386d74e61SPeter Brune 36486d74e61SPeter Brune Output Arguments: 36586d74e61SPeter Brune + Y - correction, possibly modified 36686d74e61SPeter Brune - changed - flag indicating that Y was modified 36786d74e61SPeter Brune 36886d74e61SPeter Brune Options Database Key: 369cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 370cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 37186d74e61SPeter Brune 37286d74e61SPeter Brune Level: advanced 37386d74e61SPeter Brune 37486d74e61SPeter Brune Notes: 37586d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 37686d74e61SPeter Brune 37786d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 37886d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 37986d74e61SPeter Brune is generally not useful when using a Newton linearization. 38086d74e61SPeter Brune 38186d74e61SPeter Brune Reference: 38286d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 38386d74e61SPeter Brune 38486d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 38586d74e61SPeter Brune @*/ 386f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 38786d74e61SPeter Brune { 38886d74e61SPeter Brune PetscErrorCode ierr; 38986d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 39086d74e61SPeter Brune Vec Ylast; 39186d74e61SPeter Brune PetscScalar dot; 39286d74e61SPeter Brune PetscInt iter; 39386d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 39486d74e61SPeter Brune SNES snes; 39586d74e61SPeter Brune 39686d74e61SPeter Brune PetscFunctionBegin; 397f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 39886d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 39986d74e61SPeter Brune if (!Ylast) { 40086d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 40186d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 40286d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 40386d74e61SPeter Brune } 40486d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 40586d74e61SPeter Brune if (iter < 2) { 40686d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 40786d74e61SPeter Brune *changed = PETSC_FALSE; 40886d74e61SPeter Brune PetscFunctionReturn(0); 40986d74e61SPeter Brune } 41086d74e61SPeter Brune 41186d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 41286d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 41386d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 41486d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 41586d74e61SPeter Brune theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 41686d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 41786d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 41886d74e61SPeter Brune /* Modify the step Y */ 41986d74e61SPeter Brune PetscReal alpha,ydiffnorm; 42086d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 42186d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 42286d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 42386d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 42486d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 425c69d1a72SBarry Smith ierr = PetscInfo3(snes,"Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n",(double)(theta*180./PETSC_PI),(double)angle,(double)alpha);CHKERRQ(ierr); 42686d74e61SPeter Brune } else { 427c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr); 42886d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 42986d74e61SPeter Brune *changed = PETSC_FALSE; 430bf7f4e0aSPeter Brune } 431bf7f4e0aSPeter Brune PetscFunctionReturn(0); 432bf7f4e0aSPeter Brune } 433bf7f4e0aSPeter Brune 434bf7f4e0aSPeter Brune #undef __FUNCT__ 435f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply" 436f40b411bSPeter Brune /*@ 437cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 438f40b411bSPeter Brune 439f1c6b773SPeter Brune Collective on SNESLineSearch 440f40b411bSPeter Brune 441f40b411bSPeter Brune Input Parameters: 4428e557f58SPeter Brune + linesearch - The linesearch context 4438e557f58SPeter Brune . X - The current solution 4448e557f58SPeter Brune . F - The current function 4458e557f58SPeter Brune . fnorm - The current norm 4468e557f58SPeter Brune - Y - The search direction 447f40b411bSPeter Brune 448f40b411bSPeter Brune Output Parameters: 4498e557f58SPeter Brune + X - The new solution 4508e557f58SPeter Brune . F - The new function 4518e557f58SPeter Brune - fnorm - The new function norm 452f40b411bSPeter Brune 453cd7522eaSPeter Brune Options Database Keys: 4543c7d6663SPeter Brune + -snes_linesearch_type - basic, bt, l2, cp, shell 455cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 4563c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 457cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 4583c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 4593c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 460cd7522eaSPeter Brune 461cd7522eaSPeter Brune Notes: 462cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 463cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 464cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 465cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 466cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 467cd7522eaSPeter Brune therefore may be fairly expensive. 468cd7522eaSPeter Brune 469f40b411bSPeter Brune Level: Intermediate 470f40b411bSPeter Brune 471f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 472f40b411bSPeter Brune 473cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction() 474f40b411bSPeter Brune @*/ 475*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 476*bf388a1fSBarry Smith { 477bf7f4e0aSPeter Brune PetscErrorCode ierr; 478bf7f4e0aSPeter Brune 479*bf388a1fSBarry Smith PetscFunctionBegin; 480f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 481bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 482bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 483bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 484bf7f4e0aSPeter Brune 485bf7f4e0aSPeter Brune linesearch->success = PETSC_TRUE; 486bf7f4e0aSPeter Brune 487bf7f4e0aSPeter Brune linesearch->vec_sol = X; 488bf7f4e0aSPeter Brune linesearch->vec_update = Y; 489bf7f4e0aSPeter Brune linesearch->vec_func = F; 490bf7f4e0aSPeter Brune 491f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 492bf7f4e0aSPeter Brune 493bf7f4e0aSPeter Brune if (!linesearch->keeplambda) 494bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 495bf7f4e0aSPeter Brune 496bf7f4e0aSPeter Brune if (fnorm) { 497bf7f4e0aSPeter Brune linesearch->fnorm = *fnorm; 498bf7f4e0aSPeter Brune } else { 499bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 500bf7f4e0aSPeter Brune } 501bf7f4e0aSPeter Brune 502f1c6b773SPeter Brune ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 503bf7f4e0aSPeter Brune 504bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 505bf7f4e0aSPeter Brune 506f1c6b773SPeter Brune ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 507bf7f4e0aSPeter Brune 508bf7f4e0aSPeter Brune if (fnorm) 509bf7f4e0aSPeter Brune *fnorm = linesearch->fnorm; 510bf7f4e0aSPeter Brune PetscFunctionReturn(0); 511bf7f4e0aSPeter Brune } 512bf7f4e0aSPeter Brune 513bf7f4e0aSPeter Brune #undef __FUNCT__ 514f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy" 515f40b411bSPeter Brune /*@ 516f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 517f40b411bSPeter Brune 518f1c6b773SPeter Brune Collective on SNESLineSearch 519f40b411bSPeter Brune 520f40b411bSPeter Brune Input Parameters: 5218e557f58SPeter Brune . linesearch - The linesearch context 522f40b411bSPeter Brune 523f40b411bSPeter Brune Level: Intermediate 524f40b411bSPeter Brune 52578bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 526f40b411bSPeter Brune 527cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 528f40b411bSPeter Brune @*/ 529*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 530*bf388a1fSBarry Smith { 531bf7f4e0aSPeter Brune PetscErrorCode ierr; 532*bf388a1fSBarry Smith 533bf7f4e0aSPeter Brune PetscFunctionBegin; 534bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 535f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 536bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 537bf7f4e0aSPeter Brune ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr); 53822d28d08SBarry Smith ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr); 539bf7f4e0aSPeter Brune if ((*linesearch)->ops->destroy) { 540bf7f4e0aSPeter Brune (*linesearch)->ops->destroy(*linesearch); 541bf7f4e0aSPeter Brune } 542bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 543e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 544bf7f4e0aSPeter Brune PetscFunctionReturn(0); 545bf7f4e0aSPeter Brune } 546bf7f4e0aSPeter Brune 547bf7f4e0aSPeter Brune #undef __FUNCT__ 548f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor" 549f40b411bSPeter Brune /*@ 550cd7522eaSPeter Brune SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search. 551bf7f4e0aSPeter Brune 552bf7f4e0aSPeter Brune Input Parameters: 553bf7f4e0aSPeter Brune + snes - nonlinear context obtained from SNESCreate() 554bf7f4e0aSPeter Brune - flg - PETSC_TRUE to monitor the line search 555bf7f4e0aSPeter Brune 556bf7f4e0aSPeter Brune Logically Collective on SNES 557bf7f4e0aSPeter Brune 558bf7f4e0aSPeter Brune Options Database: 5598e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 560bf7f4e0aSPeter Brune 561bf7f4e0aSPeter Brune Level: intermediate 562bf7f4e0aSPeter Brune 563bf7f4e0aSPeter Brune 564cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer 565bf7f4e0aSPeter Brune @*/ 566f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg) 567bf7f4e0aSPeter Brune { 568bf7f4e0aSPeter Brune PetscErrorCode ierr; 569*bf388a1fSBarry Smith 570bf7f4e0aSPeter Brune PetscFunctionBegin; 571bf7f4e0aSPeter Brune if (flg && !linesearch->monitor) { 572bf7f4e0aSPeter Brune ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr); 573bf7f4e0aSPeter Brune } else if (!flg && linesearch->monitor) { 574bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 575bf7f4e0aSPeter Brune } 576bf7f4e0aSPeter Brune PetscFunctionReturn(0); 577bf7f4e0aSPeter Brune } 578bf7f4e0aSPeter Brune 579bf7f4e0aSPeter Brune #undef __FUNCT__ 580f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor" 581f40b411bSPeter Brune /*@ 582cd7522eaSPeter Brune SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor. 5836a388c36SPeter Brune 584f40b411bSPeter Brune Input Parameters: 5858e557f58SPeter Brune . linesearch - linesearch context 586f40b411bSPeter Brune 587f40b411bSPeter Brune Input Parameters: 5888e557f58SPeter Brune . monitor - monitor context 589f40b411bSPeter Brune 590f40b411bSPeter Brune Logically Collective on SNES 591f40b411bSPeter Brune 592f40b411bSPeter Brune 593f40b411bSPeter Brune Options Database Keys: 5948e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 595f40b411bSPeter Brune 596f40b411bSPeter Brune Level: intermediate 597f40b411bSPeter Brune 598f40b411bSPeter Brune 599cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer 600f40b411bSPeter Brune @*/ 601f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 6026a388c36SPeter Brune { 6036a388c36SPeter Brune PetscFunctionBegin; 604f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 6056a388c36SPeter Brune if (monitor) { 6066a388c36SPeter Brune PetscValidPointer(monitor, 2); 6076a388c36SPeter Brune *monitor = linesearch->monitor; 6086a388c36SPeter Brune } 6096a388c36SPeter Brune PetscFunctionReturn(0); 6106a388c36SPeter Brune } 6116a388c36SPeter Brune 6126a388c36SPeter Brune #undef __FUNCT__ 613f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions" 614f40b411bSPeter Brune /*@ 615f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 616f40b411bSPeter Brune 617f40b411bSPeter Brune Input Parameters: 6188e557f58SPeter Brune . linesearch - linesearch context 619f40b411bSPeter Brune 620f40b411bSPeter Brune Options Database Keys: 6213c7d6663SPeter Brune + -snes_linesearch_type <type> - basic, bt, l2, cp, shell 6223c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 6235a9b6599SPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch type 62471eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 6251a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 6261a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 6271a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 6281a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 6291a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 630cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 6318e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 632cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 6331a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 6341a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method 635f40b411bSPeter Brune 636f1c6b773SPeter Brune Logically Collective on SNESLineSearch 637f40b411bSPeter Brune 638f40b411bSPeter Brune Level: intermediate 639f40b411bSPeter Brune 6403c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard() 641f40b411bSPeter Brune @*/ 642*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 643*bf388a1fSBarry Smith { 644bf7f4e0aSPeter Brune PetscErrorCode ierr; 6451a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 646bf7f4e0aSPeter Brune char type[256]; 647bf7f4e0aSPeter Brune PetscBool flg, set; 648*bf388a1fSBarry Smith 649bf7f4e0aSPeter Brune PetscFunctionBegin; 650f1c6b773SPeter Brune if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 651bf7f4e0aSPeter Brune 652bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 653bf7f4e0aSPeter Brune if (((PetscObject)linesearch)->type_name) { 654bf7f4e0aSPeter Brune deft = ((PetscObject)linesearch)->type_name; 655bf7f4e0aSPeter Brune } 6563c7d6663SPeter Brune ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 657bf7f4e0aSPeter Brune if (flg) { 658f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 659bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 660f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 661bf7f4e0aSPeter Brune } 662bf7f4e0aSPeter Brune 6637a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor", 664bf7f4e0aSPeter Brune linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 665f1c6b773SPeter Brune if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 666bf7f4e0aSPeter Brune 6671a9b3a06SPeter Brune /* tolerances */ 66871eef1aeSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr); 6691a9b3a06SPeter Brune ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr); 6701a9b3a06SPeter Brune ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr); 6711a9b3a06SPeter Brune ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr); 6721a9b3a06SPeter Brune ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr); 6738e557f58SPeter Brune ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr); 6747a35526eSPeter Brune 6751a9b3a06SPeter Brune /* damping parameters */ 6761a9b3a06SPeter Brune ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr); 6771a9b3a06SPeter Brune 6781a9b3a06SPeter Brune ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr); 6791a9b3a06SPeter Brune 6801a9b3a06SPeter Brune /* precheck */ 6817a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 6827a35526eSPeter Brune if (set) { 6837a35526eSPeter Brune if (flg) { 6847a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 6857a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 6867a35526eSPeter Brune "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr); 6877a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 6887a35526eSPeter Brune } else { 6897a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 6907a35526eSPeter Brune } 6917a35526eSPeter Brune } 692b000cd8dSPeter Brune ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr); 6931a9b3a06SPeter Brune ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr); 6947a35526eSPeter Brune 6955a9b6599SPeter Brune if (linesearch->ops->setfromoptions) { 6965a9b6599SPeter Brune (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr); 6975a9b6599SPeter Brune } 6985a9b6599SPeter Brune 699bf7f4e0aSPeter Brune ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr); 700bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 701bf7f4e0aSPeter Brune PetscFunctionReturn(0); 702bf7f4e0aSPeter Brune } 703bf7f4e0aSPeter Brune 704bf7f4e0aSPeter Brune #undef __FUNCT__ 705f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView" 706f40b411bSPeter Brune /*@ 707cd7522eaSPeter Brune SNESLineSearchView - Prints useful information about the line search not 708cd7522eaSPeter Brune related to an individual call. 709f40b411bSPeter Brune 710f40b411bSPeter Brune Input Parameters: 7118e557f58SPeter Brune . linesearch - linesearch context 712f40b411bSPeter Brune 713f1c6b773SPeter Brune Logically Collective on SNESLineSearch 714f40b411bSPeter Brune 715f40b411bSPeter Brune Level: intermediate 716f40b411bSPeter Brune 717f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 718f40b411bSPeter Brune @*/ 719*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 720*bf388a1fSBarry Smith { 7217f1410a3SPeter Brune PetscErrorCode ierr; 7227f1410a3SPeter Brune PetscBool iascii; 723*bf388a1fSBarry Smith 724bf7f4e0aSPeter Brune PetscFunctionBegin; 7257f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 7267f1410a3SPeter Brune if (!viewer) { 7277f1410a3SPeter Brune ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr); 7287f1410a3SPeter Brune } 7297f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 7307f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 731f40b411bSPeter Brune 732251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7337f1410a3SPeter Brune if (iascii) { 7347f1410a3SPeter Brune ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr); 7357f1410a3SPeter Brune if (linesearch->ops->view) { 7367f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7377f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 7387f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7397f1410a3SPeter Brune } 740c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr); 741c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr); 7427f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 7437f1410a3SPeter Brune if (linesearch->ops->precheckstep) { 7447f1410a3SPeter Brune if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) { 7457f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 7467f1410a3SPeter Brune } else { 7477f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 7487f1410a3SPeter Brune } 7497f1410a3SPeter Brune } 7507f1410a3SPeter Brune if (linesearch->ops->postcheckstep) { 7517f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 7527f1410a3SPeter Brune } 7537f1410a3SPeter Brune } 754bf7f4e0aSPeter Brune PetscFunctionReturn(0); 755bf7f4e0aSPeter Brune } 756bf7f4e0aSPeter Brune 757bf7f4e0aSPeter Brune #undef __FUNCT__ 758f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType" 759ea5d4fccSPeter Brune /*@C 760f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 761f40b411bSPeter Brune 762f40b411bSPeter Brune Input Parameters: 7638e557f58SPeter Brune + linesearch - linesearch context 764f40b411bSPeter Brune - type - The type of line search to be used 765f40b411bSPeter Brune 766cd7522eaSPeter Brune Available Types: 767cd7522eaSPeter Brune + basic - Simple damping line search. 7688e557f58SPeter Brune . bt - Backtracking line search over the L2 norm of the function 7698e557f58SPeter Brune . l2 - Secant line search over the L2 norm of the function 7708e557f58SPeter Brune . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 7718e557f58SPeter Brune - shell - User provided SNESLineSearch implementation 772cd7522eaSPeter Brune 773f1c6b773SPeter Brune Logically Collective on SNESLineSearch 774f40b411bSPeter Brune 775f40b411bSPeter Brune Level: intermediate 776f40b411bSPeter Brune 777f40b411bSPeter Brune 778f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 779f40b411bSPeter Brune @*/ 78019fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 781bf7f4e0aSPeter Brune { 782f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 783bf7f4e0aSPeter Brune PetscBool match; 784bf7f4e0aSPeter Brune 785bf7f4e0aSPeter Brune PetscFunctionBegin; 786f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 787bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 788bf7f4e0aSPeter Brune 789251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 790bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 791bf7f4e0aSPeter Brune 792f1c6b773SPeter Brune ierr = PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 793bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 794bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 795bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 796bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 797bf7f4e0aSPeter Brune linesearch->ops->destroy = PETSC_NULL; 798bf7f4e0aSPeter Brune } 799f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 800bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 801bf7f4e0aSPeter Brune linesearch->ops->view = 0; 802bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 803bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 804bf7f4e0aSPeter Brune 805bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 806bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 807bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS) 808bf7f4e0aSPeter Brune if (PetscAMSPublishAll) { 809bf7f4e0aSPeter Brune ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 810bf7f4e0aSPeter Brune } 811bf7f4e0aSPeter Brune #endif 812bf7f4e0aSPeter Brune PetscFunctionReturn(0); 813bf7f4e0aSPeter Brune } 814bf7f4e0aSPeter Brune 815bf7f4e0aSPeter Brune #undef __FUNCT__ 816f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES" 817f40b411bSPeter Brune /*@ 81878bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 819f40b411bSPeter Brune 820f40b411bSPeter Brune Input Parameters: 8218e557f58SPeter Brune + linesearch - linesearch context 822f40b411bSPeter Brune - snes - The snes instance 823f40b411bSPeter Brune 82478bcb3b5SPeter Brune Level: developer 82578bcb3b5SPeter Brune 82678bcb3b5SPeter Brune Notes: 82778bcb3b5SPeter Brune This happens automatically when the line search is gotten/created with 82878bcb3b5SPeter Brune SNESGetSNESLineSearch(). This routine is therefore mainly called within SNES 82978bcb3b5SPeter Brune implementations. 830f40b411bSPeter Brune 8318141a3b9SPeter Brune Level: developer 832f40b411bSPeter Brune 833cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 834f40b411bSPeter Brune @*/ 835*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 836*bf388a1fSBarry Smith { 837bf7f4e0aSPeter Brune PetscFunctionBegin; 838f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 839bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 840bf7f4e0aSPeter Brune linesearch->snes = snes; 841bf7f4e0aSPeter Brune PetscFunctionReturn(0); 842bf7f4e0aSPeter Brune } 843bf7f4e0aSPeter Brune 844bf7f4e0aSPeter Brune #undef __FUNCT__ 845f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES" 846f40b411bSPeter Brune /*@ 8478141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 8488141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 8498141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 8508141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 851f40b411bSPeter Brune 852f40b411bSPeter Brune Input Parameters: 8538e557f58SPeter Brune . linesearch - linesearch context 854f40b411bSPeter Brune 855f40b411bSPeter Brune Output Parameters: 856f40b411bSPeter Brune . snes - The snes instance 857f40b411bSPeter Brune 8588141a3b9SPeter Brune Level: developer 859f40b411bSPeter Brune 860cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 861f40b411bSPeter Brune @*/ 862*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 863*bf388a1fSBarry Smith { 864bf7f4e0aSPeter Brune PetscFunctionBegin; 865f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8666a388c36SPeter Brune PetscValidPointer(snes, 2); 867bf7f4e0aSPeter Brune *snes = linesearch->snes; 868bf7f4e0aSPeter Brune PetscFunctionReturn(0); 869bf7f4e0aSPeter Brune } 870bf7f4e0aSPeter Brune 8716a388c36SPeter Brune #undef __FUNCT__ 872f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda" 873f40b411bSPeter Brune /*@ 874f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 875f40b411bSPeter Brune 876f40b411bSPeter Brune Input Parameters: 8778e557f58SPeter Brune . linesearch - linesearch context 878f40b411bSPeter Brune 879f40b411bSPeter Brune Output Parameters: 880cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 881f40b411bSPeter Brune 88278bcb3b5SPeter Brune Level: advanced 88378bcb3b5SPeter Brune 8848e557f58SPeter Brune Notes: 8858e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 88678bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 88778bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 88878bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 88978bcb3b5SPeter Brune 890f40b411bSPeter Brune 891cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 892f40b411bSPeter Brune @*/ 893f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 8946a388c36SPeter Brune { 8956a388c36SPeter Brune PetscFunctionBegin; 896f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8976a388c36SPeter Brune PetscValidPointer(lambda, 2); 8986a388c36SPeter Brune *lambda = linesearch->lambda; 8996a388c36SPeter Brune PetscFunctionReturn(0); 9006a388c36SPeter Brune } 9016a388c36SPeter Brune 9026a388c36SPeter Brune #undef __FUNCT__ 903f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda" 904f40b411bSPeter Brune /*@ 905f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 906f40b411bSPeter Brune 907f40b411bSPeter Brune Input Parameters: 9088e557f58SPeter Brune + linesearch - linesearch context 909f40b411bSPeter Brune - lambda - The last steplength. 910f40b411bSPeter Brune 911cd7522eaSPeter Brune Notes: 912cd7522eaSPeter Brune This routine is typically used within implementations of SNESLineSearchApply 913cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 914cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 915cd7522eaSPeter Brune as an inner scaling parameter. 916cd7522eaSPeter Brune 91778bcb3b5SPeter Brune Level: advanced 918f40b411bSPeter Brune 919f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 920f40b411bSPeter Brune @*/ 921f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 9226a388c36SPeter Brune { 9236a388c36SPeter Brune PetscFunctionBegin; 924f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9256a388c36SPeter Brune linesearch->lambda = lambda; 9266a388c36SPeter Brune PetscFunctionReturn(0); 9276a388c36SPeter Brune } 9286a388c36SPeter Brune 9296a388c36SPeter Brune #undef __FUNCT__ 930f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances" 931f40b411bSPeter Brune /*@ 9323c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 93378bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 93478bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 93578bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 936f40b411bSPeter Brune 937f40b411bSPeter Brune Input Parameters: 9388e557f58SPeter Brune . linesearch - linesearch context 939f40b411bSPeter Brune 940f40b411bSPeter Brune Output Parameters: 941516fe3c3SPeter Brune + steptol - The minimum steplength 9426cc8e53bSPeter Brune . maxstep - The maximum steplength 943516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 944516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 945516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 946516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 947f40b411bSPeter Brune 94878bcb3b5SPeter Brune Level: intermediate 94978bcb3b5SPeter Brune 95078bcb3b5SPeter Brune Notes: 95178bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 9523c7d6663SPeter Brune the type requires. 953516fe3c3SPeter Brune 954f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 955f40b411bSPeter Brune @*/ 956f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 9576a388c36SPeter Brune { 9586a388c36SPeter Brune PetscFunctionBegin; 959f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 960516fe3c3SPeter Brune if (steptol) { 9616a388c36SPeter Brune PetscValidPointer(steptol, 2); 9626a388c36SPeter Brune *steptol = linesearch->steptol; 963516fe3c3SPeter Brune } 964516fe3c3SPeter Brune if (maxstep) { 965516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 966516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 967516fe3c3SPeter Brune } 968516fe3c3SPeter Brune if (rtol) { 969516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 970516fe3c3SPeter Brune *rtol = linesearch->rtol; 971516fe3c3SPeter Brune } 972516fe3c3SPeter Brune if (atol) { 973516fe3c3SPeter Brune PetscValidPointer(atol, 5); 974516fe3c3SPeter Brune *atol = linesearch->atol; 975516fe3c3SPeter Brune } 976516fe3c3SPeter Brune if (ltol) { 977516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 978516fe3c3SPeter Brune *ltol = linesearch->ltol; 979516fe3c3SPeter Brune } 980516fe3c3SPeter Brune if (max_its) { 981516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 982516fe3c3SPeter Brune *max_its = linesearch->max_its; 983516fe3c3SPeter Brune } 9846a388c36SPeter Brune PetscFunctionReturn(0); 9856a388c36SPeter Brune } 9866a388c36SPeter Brune 9876a388c36SPeter Brune #undef __FUNCT__ 988f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances" 989f40b411bSPeter Brune /*@ 9903c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 99178bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 99278bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 99378bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 994f40b411bSPeter Brune 995f40b411bSPeter Brune Input Parameters: 9968e557f58SPeter Brune + linesearch - linesearch context 997516fe3c3SPeter Brune . steptol - The minimum steplength 9986cc8e53bSPeter Brune . maxstep - The maximum steplength 999516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1000516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1001516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1002516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1003f40b411bSPeter Brune 100478bcb3b5SPeter Brune Notes: 10053c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 100678bcb3b5SPeter Brune place of an argument. 1007f40b411bSPeter Brune 100878bcb3b5SPeter Brune Level: intermediate 1009516fe3c3SPeter Brune 1010f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 1011f40b411bSPeter Brune @*/ 1012f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 10136a388c36SPeter Brune { 10146a388c36SPeter Brune PetscFunctionBegin; 1015f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1016d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1017d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1018d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1019d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1020d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1021d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1022d3952184SSatish Balay 1023d3952184SSatish Balay if ( steptol!= PETSC_DEFAULT) { 1024c69d1a72SBarry Smith if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 10256a388c36SPeter Brune linesearch->steptol = steptol; 1026d3952184SSatish Balay } 1027d3952184SSatish Balay 1028d3952184SSatish Balay if ( maxstep!= PETSC_DEFAULT) { 1029c69d1a72SBarry Smith if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1030516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1031d3952184SSatish Balay } 1032d3952184SSatish Balay 1033d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1034c69d1a72SBarry Smith if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol); 1035516fe3c3SPeter Brune linesearch->rtol = rtol; 1036d3952184SSatish Balay } 1037d3952184SSatish Balay 1038d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1039c69d1a72SBarry Smith if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1040516fe3c3SPeter Brune linesearch->atol = atol; 1041d3952184SSatish Balay } 1042d3952184SSatish Balay 1043d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1044c69d1a72SBarry Smith if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1045516fe3c3SPeter Brune linesearch->ltol = ltol; 1046d3952184SSatish Balay } 1047d3952184SSatish Balay 1048d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1049d3952184SSatish Balay if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1050516fe3c3SPeter Brune linesearch->max_its = max_its; 1051d3952184SSatish Balay } 1052d3952184SSatish Balay 10536a388c36SPeter Brune PetscFunctionReturn(0); 10546a388c36SPeter Brune } 10556a388c36SPeter Brune 10566a388c36SPeter Brune #undef __FUNCT__ 1057f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping" 1058f40b411bSPeter Brune /*@ 1059f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1060f40b411bSPeter Brune 1061f40b411bSPeter Brune Input Parameters: 10628e557f58SPeter Brune . linesearch - linesearch context 1063f40b411bSPeter Brune 1064f40b411bSPeter Brune Output Parameters: 10658e557f58SPeter Brune . damping - The damping parameter 1066f40b411bSPeter Brune 106778bcb3b5SPeter Brune Level: advanced 1068f40b411bSPeter Brune 106978bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1070f40b411bSPeter Brune @*/ 1071f40b411bSPeter Brune 1072f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 10736a388c36SPeter Brune { 10746a388c36SPeter Brune PetscFunctionBegin; 1075f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 10766a388c36SPeter Brune PetscValidPointer(damping, 2); 10776a388c36SPeter Brune *damping = linesearch->damping; 10786a388c36SPeter Brune PetscFunctionReturn(0); 10796a388c36SPeter Brune } 10806a388c36SPeter Brune 10816a388c36SPeter Brune #undef __FUNCT__ 1082f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping" 1083f40b411bSPeter Brune /*@ 1084f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1085f40b411bSPeter Brune 1086f40b411bSPeter Brune Input Parameters: 108778bcb3b5SPeter Brune . linesearch - linesearch context 108878bcb3b5SPeter Brune . damping - The damping parameter 1089f40b411bSPeter Brune 1090f40b411bSPeter Brune Level: intermediate 1091f40b411bSPeter Brune 1092cd7522eaSPeter Brune Notes: 1093cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1094cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 109578bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1096cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1097cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1098cd7522eaSPeter Brune 1099f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1100f40b411bSPeter Brune @*/ 1101f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 11026a388c36SPeter Brune { 11036a388c36SPeter Brune PetscFunctionBegin; 1104f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11056a388c36SPeter Brune linesearch->damping = damping; 11066a388c36SPeter Brune PetscFunctionReturn(0); 11076a388c36SPeter Brune } 11086a388c36SPeter Brune 11096a388c36SPeter Brune #undef __FUNCT__ 111059405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder" 111159405d5eSPeter Brune /*@ 111259405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 111359405d5eSPeter Brune 111459405d5eSPeter Brune Input Parameters: 111578bcb3b5SPeter Brune . linesearch - linesearch context 111659405d5eSPeter Brune 111759405d5eSPeter Brune Output Parameters: 11188e557f58SPeter Brune . order - The order 111959405d5eSPeter Brune 112078bcb3b5SPeter Brune Possible Values for order: 11213c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 11223c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 11233c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 112478bcb3b5SPeter Brune 112559405d5eSPeter Brune Level: intermediate 112659405d5eSPeter Brune 112759405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 112859405d5eSPeter Brune @*/ 112959405d5eSPeter Brune 1130b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 113159405d5eSPeter Brune { 113259405d5eSPeter Brune PetscFunctionBegin; 113359405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 113459405d5eSPeter Brune PetscValidPointer(order, 2); 113559405d5eSPeter Brune *order = linesearch->order; 113659405d5eSPeter Brune PetscFunctionReturn(0); 113759405d5eSPeter Brune } 113859405d5eSPeter Brune 113959405d5eSPeter Brune #undef __FUNCT__ 114059405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder" 114159405d5eSPeter Brune /*@ 114259405d5eSPeter Brune SNESLineSearchSetOrder - Sets the line search damping paramter. 114359405d5eSPeter Brune 114459405d5eSPeter Brune Input Parameters: 114578bcb3b5SPeter Brune . linesearch - linesearch context 114678bcb3b5SPeter Brune . order - The damping parameter 114759405d5eSPeter Brune 114859405d5eSPeter Brune Level: intermediate 114959405d5eSPeter Brune 115078bcb3b5SPeter Brune Possible Values for order: 11513c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 11523c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 11533c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 115478bcb3b5SPeter Brune 115559405d5eSPeter Brune Notes: 115659405d5eSPeter Brune Variable orders are supported by the following line searches: 115778bcb3b5SPeter Brune + bt - cubic and quadratic 115878bcb3b5SPeter Brune - cp - linear and quadratic 115959405d5eSPeter Brune 116059405d5eSPeter Brune .seealso: SNESLineSearchGetOrder() 116159405d5eSPeter Brune @*/ 1162b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 116359405d5eSPeter Brune { 116459405d5eSPeter Brune PetscFunctionBegin; 116559405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 116659405d5eSPeter Brune linesearch->order = order; 116759405d5eSPeter Brune PetscFunctionReturn(0); 116859405d5eSPeter Brune } 116959405d5eSPeter Brune 117059405d5eSPeter Brune #undef __FUNCT__ 1171f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms" 1172f40b411bSPeter Brune /*@ 1173f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1174f40b411bSPeter Brune 1175f40b411bSPeter Brune Input Parameters: 117678bcb3b5SPeter Brune . linesearch - linesearch context 1177f40b411bSPeter Brune 1178f40b411bSPeter Brune Output Parameters: 1179f40b411bSPeter Brune + xnorm - The norm of the current solution 1180f40b411bSPeter Brune . fnorm - The norm of the current function 1181f40b411bSPeter Brune - ynorm - The norm of the current update 1182f40b411bSPeter Brune 1183cd7522eaSPeter Brune Notes: 1184cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1185cd7522eaSPeter Brune 118678bcb3b5SPeter Brune Level: developer 1187f40b411bSPeter Brune 1188f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1189f40b411bSPeter Brune @*/ 1190f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1191bf7f4e0aSPeter Brune { 1192bf7f4e0aSPeter Brune PetscFunctionBegin; 1193f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1194bf7f4e0aSPeter Brune if (xnorm) { 1195bf7f4e0aSPeter Brune *xnorm = linesearch->xnorm; 1196bf7f4e0aSPeter Brune } 1197bf7f4e0aSPeter Brune if (fnorm) { 1198bf7f4e0aSPeter Brune *fnorm = linesearch->fnorm; 1199bf7f4e0aSPeter Brune } 1200bf7f4e0aSPeter Brune if (ynorm) { 1201bf7f4e0aSPeter Brune *ynorm = linesearch->ynorm; 1202bf7f4e0aSPeter Brune } 1203bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1204bf7f4e0aSPeter Brune } 1205bf7f4e0aSPeter Brune 1206e7058c64SPeter Brune #undef __FUNCT__ 1207f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms" 1208f40b411bSPeter Brune /*@ 1209f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1210f40b411bSPeter Brune 1211f40b411bSPeter Brune Input Parameters: 121278bcb3b5SPeter Brune + linesearch - linesearch context 1213f40b411bSPeter Brune . xnorm - The norm of the current solution 1214f40b411bSPeter Brune . fnorm - The norm of the current function 1215f40b411bSPeter Brune - ynorm - The norm of the current update 1216f40b411bSPeter Brune 121778bcb3b5SPeter Brune Level: advanced 1218f40b411bSPeter Brune 1219f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1220f40b411bSPeter Brune @*/ 1221f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 12226a388c36SPeter Brune { 12236a388c36SPeter Brune PetscFunctionBegin; 1224f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12256a388c36SPeter Brune linesearch->xnorm = xnorm; 12266a388c36SPeter Brune linesearch->fnorm = fnorm; 12276a388c36SPeter Brune linesearch->ynorm = ynorm; 12286a388c36SPeter Brune PetscFunctionReturn(0); 12296a388c36SPeter Brune } 12306a388c36SPeter Brune 12316a388c36SPeter Brune #undef __FUNCT__ 1232f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms" 1233f40b411bSPeter Brune /*@ 1234f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1235f40b411bSPeter Brune 1236f40b411bSPeter Brune Input Parameters: 123778bcb3b5SPeter Brune . linesearch - linesearch context 1238f40b411bSPeter Brune 1239f40b411bSPeter Brune Options Database Keys: 12408e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1241f40b411bSPeter Brune 1242f40b411bSPeter Brune Level: intermediate 1243f40b411bSPeter Brune 124478bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1245f40b411bSPeter Brune @*/ 1246f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 12476a388c36SPeter Brune { 12486a388c36SPeter Brune PetscErrorCode ierr; 12499bd66eb0SPeter Brune SNES snes; 1250*bf388a1fSBarry Smith 12516a388c36SPeter Brune PetscFunctionBegin; 12526a388c36SPeter Brune if (linesearch->norms) { 12539bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1254f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 12559bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12569bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12579bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 12589bd66eb0SPeter Brune } else { 12596a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12606a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12616a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12626a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12636a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12646a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12656a388c36SPeter Brune } 12669bd66eb0SPeter Brune } 12676a388c36SPeter Brune PetscFunctionReturn(0); 12686a388c36SPeter Brune } 12696a388c36SPeter Brune 12706f263ca3SPeter Brune #undef __FUNCT__ 12716f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms" 12726f263ca3SPeter Brune /*@ 12736f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 12746f263ca3SPeter Brune 12756f263ca3SPeter Brune Input Parameters: 127678bcb3b5SPeter Brune + linesearch - linesearch context 127778bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 12786f263ca3SPeter Brune 12796f263ca3SPeter Brune Options Database Keys: 12808e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 12816f263ca3SPeter Brune 12826f263ca3SPeter Brune Notes: 12831a4f838cSPeter Brune This is most relevant to the SNESLINESEARCHBASIC line search type. 12846f263ca3SPeter Brune 12856f263ca3SPeter Brune Level: intermediate 12866f263ca3SPeter Brune 12871a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 12886f263ca3SPeter Brune @*/ 12896f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 12906f263ca3SPeter Brune { 12916f263ca3SPeter Brune PetscFunctionBegin; 12926f263ca3SPeter Brune linesearch->norms = flg; 12936f263ca3SPeter Brune PetscFunctionReturn(0); 12946f263ca3SPeter Brune } 12956f263ca3SPeter Brune 12966a388c36SPeter Brune #undef __FUNCT__ 1297f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs" 1298f40b411bSPeter Brune /*@ 1299f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1300f40b411bSPeter Brune 1301f40b411bSPeter Brune Input Parameters: 130278bcb3b5SPeter Brune . linesearch - linesearch context 1303f40b411bSPeter Brune 1304f40b411bSPeter Brune Output Parameters: 1305f40b411bSPeter Brune + X - The old solution 1306f40b411bSPeter Brune . F - The old function 1307f40b411bSPeter Brune . Y - The search direction 1308f40b411bSPeter Brune . W - The new solution 1309f40b411bSPeter Brune - G - The new function 1310f40b411bSPeter Brune 131178bcb3b5SPeter Brune Level: advanced 1312f40b411bSPeter Brune 1313f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1314f40b411bSPeter Brune @*/ 1315*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1316*bf388a1fSBarry Smith { 13176a388c36SPeter Brune PetscFunctionBegin; 1318f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13196a388c36SPeter Brune if (X) { 13206a388c36SPeter Brune PetscValidPointer(X, 2); 13216a388c36SPeter Brune *X = linesearch->vec_sol; 13226a388c36SPeter Brune } 13236a388c36SPeter Brune if (F) { 13246a388c36SPeter Brune PetscValidPointer(F, 3); 13256a388c36SPeter Brune *F = linesearch->vec_func; 13266a388c36SPeter Brune } 13276a388c36SPeter Brune if (Y) { 13286a388c36SPeter Brune PetscValidPointer(Y, 4); 13296a388c36SPeter Brune *Y = linesearch->vec_update; 13306a388c36SPeter Brune } 13316a388c36SPeter Brune if (W) { 13326a388c36SPeter Brune PetscValidPointer(W, 5); 13336a388c36SPeter Brune *W = linesearch->vec_sol_new; 13346a388c36SPeter Brune } 13356a388c36SPeter Brune if (G) { 13366a388c36SPeter Brune PetscValidPointer(G, 6); 13376a388c36SPeter Brune *G = linesearch->vec_func_new; 13386a388c36SPeter Brune } 13396a388c36SPeter Brune 13406a388c36SPeter Brune PetscFunctionReturn(0); 13416a388c36SPeter Brune } 13426a388c36SPeter Brune 13436a388c36SPeter Brune #undef __FUNCT__ 1344f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs" 1345f40b411bSPeter Brune /*@ 1346f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1347f40b411bSPeter Brune 1348f40b411bSPeter Brune Input Parameters: 134978bcb3b5SPeter Brune + linesearch - linesearch context 1350f40b411bSPeter Brune . X - The old solution 1351f40b411bSPeter Brune . F - The old function 1352f40b411bSPeter Brune . Y - The search direction 1353f40b411bSPeter Brune . W - The new solution 1354f40b411bSPeter Brune - G - The new function 1355f40b411bSPeter Brune 135678bcb3b5SPeter Brune Level: advanced 1357f40b411bSPeter Brune 1358f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1359f40b411bSPeter Brune @*/ 1360*bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1361*bf388a1fSBarry Smith { 13626a388c36SPeter Brune PetscFunctionBegin; 1363f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13646a388c36SPeter Brune if (X) { 13656a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 13666a388c36SPeter Brune linesearch->vec_sol = X; 13676a388c36SPeter Brune } 13686a388c36SPeter Brune if (F) { 13696a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 13706a388c36SPeter Brune linesearch->vec_func = F; 13716a388c36SPeter Brune } 13726a388c36SPeter Brune if (Y) { 13736a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 13746a388c36SPeter Brune linesearch->vec_update = Y; 13756a388c36SPeter Brune } 13766a388c36SPeter Brune if (W) { 13776a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 13786a388c36SPeter Brune linesearch->vec_sol_new = W; 13796a388c36SPeter Brune } 13806a388c36SPeter Brune if (G) { 13816a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 13826a388c36SPeter Brune linesearch->vec_func_new = G; 13836a388c36SPeter Brune } 13846a388c36SPeter Brune PetscFunctionReturn(0); 13856a388c36SPeter Brune } 13866a388c36SPeter Brune 13876a388c36SPeter Brune #undef __FUNCT__ 1388f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1389e7058c64SPeter Brune /*@C 1390f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1391e7058c64SPeter Brune SNES options in the database. 1392e7058c64SPeter Brune 1393cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1394e7058c64SPeter Brune 1395e7058c64SPeter Brune Input Parameters: 1396e7058c64SPeter Brune + snes - the SNES context 1397e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1398e7058c64SPeter Brune 1399e7058c64SPeter Brune Notes: 1400e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1401e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1402e7058c64SPeter Brune 1403e7058c64SPeter Brune Level: advanced 1404e7058c64SPeter Brune 1405f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1406e7058c64SPeter Brune 1407e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1408e7058c64SPeter Brune @*/ 1409f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1410e7058c64SPeter Brune { 1411e7058c64SPeter Brune PetscErrorCode ierr; 1412e7058c64SPeter Brune 1413e7058c64SPeter Brune PetscFunctionBegin; 1414f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1415e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1416e7058c64SPeter Brune PetscFunctionReturn(0); 1417e7058c64SPeter Brune } 1418e7058c64SPeter Brune 1419e7058c64SPeter Brune #undef __FUNCT__ 1420f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1421e7058c64SPeter Brune /*@C 1422f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1423f1c6b773SPeter Brune SNESLineSearch options in the database. 1424e7058c64SPeter Brune 1425e7058c64SPeter Brune Not Collective 1426e7058c64SPeter Brune 1427e7058c64SPeter Brune Input Parameter: 1428cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1429e7058c64SPeter Brune 1430e7058c64SPeter Brune Output Parameter: 1431e7058c64SPeter Brune . prefix - pointer to the prefix string used 1432e7058c64SPeter Brune 14338e557f58SPeter Brune Notes: 14348e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1435e7058c64SPeter Brune sufficient length to hold the prefix. 1436e7058c64SPeter Brune 1437e7058c64SPeter Brune Level: advanced 1438e7058c64SPeter Brune 1439f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1440e7058c64SPeter Brune 1441e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1442e7058c64SPeter Brune @*/ 1443f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1444e7058c64SPeter Brune { 1445e7058c64SPeter Brune PetscErrorCode ierr; 1446e7058c64SPeter Brune 1447e7058c64SPeter Brune PetscFunctionBegin; 1448f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1449e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1450e7058c64SPeter Brune PetscFunctionReturn(0); 1451e7058c64SPeter Brune } 1452bf7f4e0aSPeter Brune 1453bf7f4e0aSPeter Brune #undef __FUNCT__ 1454f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork" 1455f40b411bSPeter Brune /*@ 1456f1c6b773SPeter Brune SNESLineSearchGetWork - Gets work vectors for the line search. 1457f40b411bSPeter Brune 1458f40b411bSPeter Brune Input Parameter: 1459f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1460f40b411bSPeter Brune - nwork - the number of work vectors 1461f40b411bSPeter Brune 1462f40b411bSPeter Brune Level: developer 1463f40b411bSPeter Brune 1464cd7522eaSPeter Brune Notes: 1465cd7522eaSPeter Brune This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation. 1466cd7522eaSPeter Brune 1467f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1468f40b411bSPeter Brune 1469f40b411bSPeter Brune .seealso: SNESDefaultGetWork() 1470f40b411bSPeter Brune @*/ 1471f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork) 1472bf7f4e0aSPeter Brune { 1473bf7f4e0aSPeter Brune PetscErrorCode ierr; 1474*bf388a1fSBarry Smith 1475bf7f4e0aSPeter Brune PetscFunctionBegin; 1476bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1477bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1478*bf388a1fSBarry Smith } SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1479bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1480bf7f4e0aSPeter Brune } 1481bf7f4e0aSPeter Brune 1482bf7f4e0aSPeter Brune #undef __FUNCT__ 1483f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess" 1484f40b411bSPeter Brune /*@ 1485f1c6b773SPeter Brune SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application 1486f40b411bSPeter Brune 1487f40b411bSPeter Brune Input Parameters: 148878bcb3b5SPeter Brune . linesearch - linesearch context 1489f40b411bSPeter Brune 1490f40b411bSPeter Brune Output Parameters: 14918e557f58SPeter Brune . success - The success or failure status 1492f40b411bSPeter Brune 1493cd7522eaSPeter Brune Notes: 1494cd7522eaSPeter Brune This is typically called after SNESLineSearchApply in order to determine if the line-search failed 1495cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1496cd7522eaSPeter Brune 1497f40b411bSPeter Brune Level: intermediate 1498f40b411bSPeter Brune 1499f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess() 1500f40b411bSPeter Brune @*/ 1501f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success) 1502bf7f4e0aSPeter Brune { 1503bf7f4e0aSPeter Brune PetscFunctionBegin; 1504f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15056a388c36SPeter Brune PetscValidPointer(success, 2); 1506bf7f4e0aSPeter Brune if (success) { 1507bf7f4e0aSPeter Brune *success = linesearch->success; 1508bf7f4e0aSPeter Brune } 1509bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1510bf7f4e0aSPeter Brune } 1511bf7f4e0aSPeter Brune 1512bf7f4e0aSPeter Brune #undef __FUNCT__ 1513f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess" 1514f40b411bSPeter Brune /*@ 1515f1c6b773SPeter Brune SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application 1516f40b411bSPeter Brune 1517f40b411bSPeter Brune Input Parameters: 151878bcb3b5SPeter Brune + linesearch - linesearch context 15198e557f58SPeter Brune - success - The success or failure status 1520f40b411bSPeter Brune 1521cd7522eaSPeter Brune Notes: 1522cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1523cd7522eaSPeter Brune the success or failure of the line search method. 1524cd7522eaSPeter Brune 152578bcb3b5SPeter Brune Level: developer 1526f40b411bSPeter Brune 1527f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess() 1528f40b411bSPeter Brune @*/ 1529f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success) 15306a388c36SPeter Brune { 1531f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15326a388c36SPeter Brune PetscFunctionBegin; 15336a388c36SPeter Brune linesearch->success = success; 15346a388c36SPeter Brune PetscFunctionReturn(0); 15356a388c36SPeter Brune } 15366a388c36SPeter Brune 15376a388c36SPeter Brune #undef __FUNCT__ 1538f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions" 15399bd66eb0SPeter Brune /*@C 1540f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 15419bd66eb0SPeter Brune 15429bd66eb0SPeter Brune Input Parameters: 15439bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 15449bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 15459bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 15469bd66eb0SPeter Brune 15479bd66eb0SPeter Brune Logically Collective on SNES 15489bd66eb0SPeter Brune 15499bd66eb0SPeter Brune Calling sequence of projectfunc: 15509bd66eb0SPeter Brune .vb 15519bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 15529bd66eb0SPeter Brune .ve 15539bd66eb0SPeter Brune 15549bd66eb0SPeter Brune Input parameters for projectfunc: 15559bd66eb0SPeter Brune + snes - nonlinear context 15569bd66eb0SPeter Brune - X - current solution 15579bd66eb0SPeter Brune 1558cd7522eaSPeter Brune Output parameters for projectfunc: 15599bd66eb0SPeter Brune . X - Projected solution 15609bd66eb0SPeter Brune 15619bd66eb0SPeter Brune Calling sequence of normfunc: 15629bd66eb0SPeter Brune .vb 15639bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 15649bd66eb0SPeter Brune .ve 15659bd66eb0SPeter Brune 1566cd7522eaSPeter Brune Input parameters for normfunc: 15679bd66eb0SPeter Brune + snes - nonlinear context 15689bd66eb0SPeter Brune . X - current solution 15699bd66eb0SPeter Brune - F - current residual 15709bd66eb0SPeter Brune 1571cd7522eaSPeter Brune Output parameters for normfunc: 15729bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 15739bd66eb0SPeter Brune 1574cd7522eaSPeter Brune Notes: 1575cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1576cd7522eaSPeter Brune 1577cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1578cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 15799bd66eb0SPeter Brune 15809bd66eb0SPeter Brune Level: developer 15819bd66eb0SPeter Brune 15829bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 15839bd66eb0SPeter Brune 1584f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 15859bd66eb0SPeter Brune @*/ 1586f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 15879bd66eb0SPeter Brune { 15889bd66eb0SPeter Brune PetscFunctionBegin; 1589f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15909bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 15919bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 15929bd66eb0SPeter Brune PetscFunctionReturn(0); 15939bd66eb0SPeter Brune } 15949bd66eb0SPeter Brune 15959bd66eb0SPeter Brune #undef __FUNCT__ 1596f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions" 15979bd66eb0SPeter Brune /*@C 1598f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 15999bd66eb0SPeter Brune 16009bd66eb0SPeter Brune Input Parameters: 16019bd66eb0SPeter Brune . snes - nonlinear context obtained from SNESCreate() 16029bd66eb0SPeter Brune 16039bd66eb0SPeter Brune Output Parameters: 16049bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16059bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16069bd66eb0SPeter Brune 16079bd66eb0SPeter Brune Logically Collective on SNES 16089bd66eb0SPeter Brune 16099bd66eb0SPeter Brune Level: developer 16109bd66eb0SPeter Brune 16119bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 16129bd66eb0SPeter Brune 1613f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 16149bd66eb0SPeter Brune @*/ 1615f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 16169bd66eb0SPeter Brune { 16179bd66eb0SPeter Brune PetscFunctionBegin; 16189bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16199bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16209bd66eb0SPeter Brune PetscFunctionReturn(0); 16219bd66eb0SPeter Brune } 16229bd66eb0SPeter Brune 16239bd66eb0SPeter Brune #undef __FUNCT__ 1624f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister" 1625bf7f4e0aSPeter Brune /*@C 1626f1c6b773SPeter Brune SNESLineSearchRegister - See SNESLineSearchRegisterDynamic() 1627bf7f4e0aSPeter Brune 1628bf7f4e0aSPeter Brune Level: advanced 1629bf7f4e0aSPeter Brune @*/ 1630f1c6b773SPeter Brune PetscErrorCode SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch)) 1631bf7f4e0aSPeter Brune { 1632bf7f4e0aSPeter Brune char fullname[PETSC_MAX_PATH_LEN]; 1633bf7f4e0aSPeter Brune PetscErrorCode ierr; 1634bf7f4e0aSPeter Brune 1635bf7f4e0aSPeter Brune PetscFunctionBegin; 1636bf7f4e0aSPeter Brune ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1637f1c6b773SPeter Brune ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1638bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1639bf7f4e0aSPeter Brune } 1640