1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2bf7f4e0aSPeter Brune 3f1c6b773SPeter Brune PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 40298fd71SBarry Smith PetscFunctionList SNESLineSearchList = NULL; 5bf7f4e0aSPeter Brune 6f1c6b773SPeter Brune PetscClassId SNESLINESEARCH_CLASSID; 7f1c6b773SPeter Brune PetscLogEvent SNESLineSearch_Apply; 8bf7f4e0aSPeter Brune 9bf7f4e0aSPeter Brune #undef __FUNCT__ 10dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorCancel" 11dcf2fd19SBarry Smith /*@ 12dcf2fd19SBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object. 13dcf2fd19SBarry Smith 14dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 15dcf2fd19SBarry Smith 16dcf2fd19SBarry Smith Input Parameters: 17dcf2fd19SBarry Smith . ls - the SNESLineSearch context 18dcf2fd19SBarry Smith 19dcf2fd19SBarry Smith Options Database Key: 20dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 21dcf2fd19SBarry Smith into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those 22dcf2fd19SBarry Smith set via the options database 23dcf2fd19SBarry Smith 24dcf2fd19SBarry Smith Notes: 25dcf2fd19SBarry Smith There is no way to clear one specific monitor from a SNESLineSearch object. 26dcf2fd19SBarry Smith 27dcf2fd19SBarry Smith This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel 28dcf2fd19SBarry Smith that one. 29dcf2fd19SBarry Smith 30dcf2fd19SBarry Smith Level: intermediate 31dcf2fd19SBarry Smith 32dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 33dcf2fd19SBarry Smith 34dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet() 35dcf2fd19SBarry Smith @*/ 36dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 37dcf2fd19SBarry Smith { 38dcf2fd19SBarry Smith PetscErrorCode ierr; 39dcf2fd19SBarry Smith PetscInt i; 40dcf2fd19SBarry Smith 41dcf2fd19SBarry Smith PetscFunctionBegin; 42dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 43dcf2fd19SBarry Smith for (i=0; i<ls->numbermonitors; i++) { 44dcf2fd19SBarry Smith if (ls->monitordestroy[i]) { 45dcf2fd19SBarry Smith ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr); 46dcf2fd19SBarry Smith } 47dcf2fd19SBarry Smith } 48dcf2fd19SBarry Smith ls->numbermonitors = 0; 49dcf2fd19SBarry Smith PetscFunctionReturn(0); 50dcf2fd19SBarry Smith } 51dcf2fd19SBarry Smith 52dcf2fd19SBarry Smith #undef __FUNCT__ 53dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitor" 54dcf2fd19SBarry Smith /*@ 55dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 56dcf2fd19SBarry Smith 57dcf2fd19SBarry Smith Collective on SNES 58dcf2fd19SBarry Smith 59dcf2fd19SBarry Smith Input Parameters: 60dcf2fd19SBarry Smith . ls - the linesearch object 61dcf2fd19SBarry Smith 62dcf2fd19SBarry Smith Notes: 63dcf2fd19SBarry Smith This routine is called by the SNES implementations. 64dcf2fd19SBarry Smith It does not typically need to be called by the user. 65dcf2fd19SBarry Smith 66dcf2fd19SBarry Smith Level: developer 67dcf2fd19SBarry Smith 68dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorSet() 69dcf2fd19SBarry Smith @*/ 70dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 71dcf2fd19SBarry Smith { 72dcf2fd19SBarry Smith PetscErrorCode ierr; 73dcf2fd19SBarry Smith PetscInt i,n = ls->numbermonitors; 74dcf2fd19SBarry Smith 75dcf2fd19SBarry Smith PetscFunctionBegin; 76dcf2fd19SBarry Smith for (i=0; i<n; i++) { 77dcf2fd19SBarry Smith ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr); 78dcf2fd19SBarry Smith } 79dcf2fd19SBarry Smith PetscFunctionReturn(0); 80dcf2fd19SBarry Smith } 81dcf2fd19SBarry Smith 82dcf2fd19SBarry Smith #undef __FUNCT__ 83dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSet" 84dcf2fd19SBarry Smith /*@C 85dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 86dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 87dcf2fd19SBarry Smith progress. 88dcf2fd19SBarry Smith 89dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 90dcf2fd19SBarry Smith 91dcf2fd19SBarry Smith Input Parameters: 92dcf2fd19SBarry Smith + ls - the SNESLineSearch context 93dcf2fd19SBarry Smith . f - the monitor function 94dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 95dcf2fd19SBarry Smith monitor routine (use NULL if no context is desired) 96dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 97dcf2fd19SBarry Smith (may be NULL) 98dcf2fd19SBarry Smith 99dcf2fd19SBarry Smith Notes: 100dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 101dcf2fd19SBarry Smith SNESLineSearchMonitorSet() multiple times; all will be called in the 102dcf2fd19SBarry Smith order in which they were set. 103dcf2fd19SBarry Smith 104dcf2fd19SBarry Smith Fortran notes: Only a single monitor function can be set for each SNESLineSearch object 105dcf2fd19SBarry Smith 106dcf2fd19SBarry Smith Level: intermediate 107dcf2fd19SBarry Smith 108dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 109dcf2fd19SBarry Smith 110dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel() 111dcf2fd19SBarry Smith @*/ 112dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 113dcf2fd19SBarry Smith { 114*78064530SBarry Smith PetscErrorCode ierr; 115*78064530SBarry Smith PetscInt i; 116*78064530SBarry Smith PetscBool identical; 117*78064530SBarry Smith 118dcf2fd19SBarry Smith PetscFunctionBegin; 119dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 120*78064530SBarry Smith for (i=0; i<ls->numbermonitors;i++) { 121*78064530SBarry Smith ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr); 122*78064530SBarry Smith if (identical) PetscFunctionReturn(0); 123*78064530SBarry Smith } 124dcf2fd19SBarry Smith if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 125dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 126dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 127dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void*)mctx; 128dcf2fd19SBarry Smith PetscFunctionReturn(0); 129dcf2fd19SBarry Smith } 130dcf2fd19SBarry Smith 131dcf2fd19SBarry Smith #undef __FUNCT__ 132dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSolutionUpdate" 133dcf2fd19SBarry Smith /*@C 134dcf2fd19SBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries 135dcf2fd19SBarry Smith 136dcf2fd19SBarry Smith Collective on SNESLineSearch 137dcf2fd19SBarry Smith 138dcf2fd19SBarry Smith Input Parameters: 139dcf2fd19SBarry Smith + ls - the SNES linesearch object 140d12e167eSBarry Smith - vf - the context for the monitor, in this case it is an ASCII PetscViewer and format 141dcf2fd19SBarry Smith 142dcf2fd19SBarry Smith Level: intermediate 143dcf2fd19SBarry Smith 144dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm 145dcf2fd19SBarry Smith 146dcf2fd19SBarry Smith .seealso: SNESMonitorSet(), SNESMonitorSolution() 147dcf2fd19SBarry Smith @*/ 148d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf) 149dcf2fd19SBarry Smith { 150dcf2fd19SBarry Smith PetscErrorCode ierr; 151d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 152dcf2fd19SBarry Smith Vec Y,W,G; 153dcf2fd19SBarry Smith 154dcf2fd19SBarry Smith PetscFunctionBegin; 155dcf2fd19SBarry Smith ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr); 156d12e167eSBarry Smith ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 157dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr); 158dcf2fd19SBarry Smith ierr = VecView(Y,viewer);CHKERRQ(ierr); 159dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr); 160dcf2fd19SBarry Smith ierr = VecView(W,viewer);CHKERRQ(ierr); 161dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr); 162dcf2fd19SBarry Smith ierr = VecView(G,viewer);CHKERRQ(ierr); 163d12e167eSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 164dcf2fd19SBarry Smith PetscFunctionReturn(0); 165dcf2fd19SBarry Smith } 166dcf2fd19SBarry Smith 167dcf2fd19SBarry Smith #undef __FUNCT__ 168f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate" 169f40b411bSPeter Brune /*@ 170cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 171f40b411bSPeter Brune 172cd7522eaSPeter Brune Logically Collective on Comm 173f40b411bSPeter Brune 174f40b411bSPeter Brune Input Parameters: 175cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 176f40b411bSPeter Brune 177f40b411bSPeter Brune Output Parameters: 1788e557f58SPeter Brune . outlinesearch - the new linesearch context 179f40b411bSPeter Brune 180162e0bf5SPeter Brune Level: developer 181162e0bf5SPeter Brune 182162e0bf5SPeter Brune Notes: 183162e0bf5SPeter Brune The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance 184162e0bf5SPeter Brune already associated with the SNES. This function is for developer use. 185f40b411bSPeter Brune 186cd7522eaSPeter Brune .keywords: LineSearch, create, context 187f40b411bSPeter Brune 188162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch() 189f40b411bSPeter Brune @*/ 190f40b411bSPeter Brune 191bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 192bf388a1fSBarry Smith { 193bf7f4e0aSPeter Brune PetscErrorCode ierr; 194f1c6b773SPeter Brune SNESLineSearch linesearch; 195bf388a1fSBarry Smith 196bf7f4e0aSPeter Brune PetscFunctionBegin; 197ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 198f34a81feSMatthew G. Knepley ierr = SNESInitializePackage();CHKERRQ(ierr); 1990298fd71SBarry Smith *outlinesearch = NULL; 200f5af7f23SKarl Rupp 20173107ff1SLisandro Dalcin ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 202bf7f4e0aSPeter Brune 2030298fd71SBarry Smith linesearch->vec_sol_new = NULL; 2040298fd71SBarry Smith linesearch->vec_func_new = NULL; 2050298fd71SBarry Smith linesearch->vec_sol = NULL; 2060298fd71SBarry Smith linesearch->vec_func = NULL; 2070298fd71SBarry Smith linesearch->vec_update = NULL; 2089bd66eb0SPeter Brune 209bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 210bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 211bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 212bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 213422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 214bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 215bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 216bf7f4e0aSPeter Brune linesearch->damping = 1.0; 217bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 218bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 219516fe3c3SPeter Brune linesearch->rtol = 1e-8; 220516fe3c3SPeter Brune linesearch->atol = 1e-15; 221516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2220298fd71SBarry Smith linesearch->precheckctx = NULL; 2230298fd71SBarry Smith linesearch->postcheckctx = NULL; 22459405d5eSPeter Brune linesearch->max_its = 1; 225bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 226bf7f4e0aSPeter Brune *outlinesearch = linesearch; 227bf7f4e0aSPeter Brune PetscFunctionReturn(0); 228bf7f4e0aSPeter Brune } 229bf7f4e0aSPeter Brune 230bf7f4e0aSPeter Brune #undef __FUNCT__ 231f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp" 232f40b411bSPeter Brune /*@ 23378bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 23478bcb3b5SPeter Brune any required vectors. 235f40b411bSPeter Brune 236cd7522eaSPeter Brune Collective on SNESLineSearch 237f40b411bSPeter Brune 238f40b411bSPeter Brune Input Parameters: 239f40b411bSPeter Brune . linesearch - The LineSearch instance. 240f40b411bSPeter Brune 241cd7522eaSPeter Brune Notes: 242f190f2fcSBarry Smith For most cases, this needn't be called by users or outside of SNESLineSearchApply(). 243cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 24478bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 245cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 246cd7522eaSPeter Brune allocated upfront. 247cd7522eaSPeter Brune 24878bcb3b5SPeter Brune Level: advanced 249f40b411bSPeter Brune 250f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 251f40b411bSPeter Brune 252f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 253f40b411bSPeter Brune @*/ 254f40b411bSPeter Brune 255bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 256bf388a1fSBarry Smith { 257bf7f4e0aSPeter Brune PetscErrorCode ierr; 258bf388a1fSBarry Smith 259bf7f4e0aSPeter Brune PetscFunctionBegin; 260bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 2611a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 262bf7f4e0aSPeter Brune } 263bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 2649bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 265bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 2669bd66eb0SPeter Brune } 2679bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 2689bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 2699bd66eb0SPeter Brune } 270bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 271bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 272bf7f4e0aSPeter Brune } 273ed07d7d7SPeter Brune if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);} 274bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 275bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 276bf7f4e0aSPeter Brune } 277bf7f4e0aSPeter Brune PetscFunctionReturn(0); 278bf7f4e0aSPeter Brune } 279bf7f4e0aSPeter Brune 280bf7f4e0aSPeter Brune #undef __FUNCT__ 281f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset" 282f40b411bSPeter Brune 283f40b411bSPeter Brune /*@ 284f190f2fcSBarry Smith SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search. 285f40b411bSPeter Brune 286f1c6b773SPeter Brune Collective on SNESLineSearch 287f40b411bSPeter Brune 288f40b411bSPeter Brune Input Parameters: 289f40b411bSPeter Brune . linesearch - The LineSearch instance. 290f40b411bSPeter Brune 291f190f2fcSBarry Smith Notes: Usually only called by SNESReset() 292f190f2fcSBarry Smith 293f190f2fcSBarry Smith Level: developer 294f40b411bSPeter Brune 295cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 296f40b411bSPeter Brune 297f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 298f40b411bSPeter Brune @*/ 299f40b411bSPeter Brune 300bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 301bf388a1fSBarry Smith { 302bf7f4e0aSPeter Brune PetscErrorCode ierr; 303bf388a1fSBarry Smith 304bf7f4e0aSPeter Brune PetscFunctionBegin; 305f5af7f23SKarl Rupp if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch); 306f5af7f23SKarl Rupp 307bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 308bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 309bf7f4e0aSPeter Brune 310bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 311f5af7f23SKarl Rupp 312bf7f4e0aSPeter Brune linesearch->nwork = 0; 313bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 314bf7f4e0aSPeter Brune PetscFunctionReturn(0); 315bf7f4e0aSPeter Brune } 316bf7f4e0aSPeter Brune 317ed07d7d7SPeter Brune #undef __FUNCT__ 318ed07d7d7SPeter Brune #define __FUNCT__ "SNESLineSearchSetFunction" 319ed07d7d7SPeter Brune /*@C 320f190f2fcSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search 321ed07d7d7SPeter Brune 322ed07d7d7SPeter Brune Input Parameters: 323ed07d7d7SPeter Brune . linesearch - the SNESLineSearch context 324f190f2fcSBarry Smith + func - function evaluation routine 325ed07d7d7SPeter Brune 326ed07d7d7SPeter Brune Level: developer 327ed07d7d7SPeter Brune 328f190f2fcSBarry Smith Notes: This is used internally by PETSc and not called by users 329f190f2fcSBarry Smith 330ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check 331ed07d7d7SPeter Brune 332f190f2fcSBarry Smith .seealso: SNESSetFunction() 333ed07d7d7SPeter Brune @*/ 334ed07d7d7SPeter Brune PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec)) 335ed07d7d7SPeter Brune { 336ed07d7d7SPeter Brune PetscFunctionBegin; 337ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 338ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 339ed07d7d7SPeter Brune PetscFunctionReturn(0); 340ed07d7d7SPeter Brune } 341ed07d7d7SPeter Brune 342ed07d7d7SPeter Brune 3436b2b7091SBarry Smith /*MC 344f190f2fcSBarry Smith SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called 3456b2b7091SBarry Smith 3466b2b7091SBarry Smith Synopsis: 347aaa7dc30SBarry Smith #include <petscsnes.h> 3486b2b7091SBarry Smith SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 3496b2b7091SBarry Smith 3506b2b7091SBarry Smith Input Parameters: 3516b2b7091SBarry Smith + x - solution vector 3526b2b7091SBarry Smith . y - search direction vector 3536b2b7091SBarry Smith - changed - flag to indicate the precheck changed x or y. 3546b2b7091SBarry Smith 355f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck() 356f190f2fcSBarry Smith and SNESLineSearchGetPreCheck() 357f190f2fcSBarry Smith 358878cb397SSatish Balay Level: advanced 359878cb397SSatish Balay 360f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck() 3616b2b7091SBarry Smith M*/ 3626b2b7091SBarry Smith 36386d74e61SPeter Brune #undef __FUNCT__ 364f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck" 36586d74e61SPeter Brune /*@C 366f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 367df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 368f190f2fcSBarry Smith determined the search direction. 36986d74e61SPeter Brune 370f1c6b773SPeter Brune Logically Collective on SNESLineSearch 37186d74e61SPeter Brune 37286d74e61SPeter Brune Input Parameters: 373f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 374f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence 375f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 37686d74e61SPeter Brune 37786d74e61SPeter Brune Level: intermediate 37886d74e61SPeter Brune 37986d74e61SPeter Brune .keywords: set, linesearch, pre-check 38086d74e61SPeter Brune 381f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 38286d74e61SPeter Brune @*/ 383f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx) 38486d74e61SPeter Brune { 3859bd66eb0SPeter Brune PetscFunctionBegin; 386f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 387f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 38886d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 38986d74e61SPeter Brune PetscFunctionReturn(0); 39086d74e61SPeter Brune } 39186d74e61SPeter Brune 39286d74e61SPeter Brune #undef __FUNCT__ 393f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck" 39486d74e61SPeter Brune /*@C 39553deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 39686d74e61SPeter Brune 39786d74e61SPeter Brune Input Parameters: 398f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 39986d74e61SPeter Brune 40086d74e61SPeter Brune Output Parameters: 401f190f2fcSBarry Smith + func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence 402f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 40386d74e61SPeter Brune 40486d74e61SPeter Brune Level: intermediate 40586d74e61SPeter Brune 40686d74e61SPeter Brune .keywords: get, linesearch, pre-check 40786d74e61SPeter Brune 408f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 40986d74e61SPeter Brune @*/ 410f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx) 41186d74e61SPeter Brune { 4129bd66eb0SPeter Brune PetscFunctionBegin; 413f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 414f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 41586d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 41686d74e61SPeter Brune PetscFunctionReturn(0); 41786d74e61SPeter Brune } 41886d74e61SPeter Brune 4196b2b7091SBarry Smith /*MC 420f190f2fcSBarry Smith SNESLineSearchPostCheckFunction - form of function that is called after line search is complete 4216b2b7091SBarry Smith 4226b2b7091SBarry Smith Synopsis: 423aaa7dc30SBarry Smith #include <petscsnes.h> 4246b2b7091SBarry Smith SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y, Vec w, *changed_y, PetscBool *changed_w); 4256b2b7091SBarry Smith 4266b2b7091SBarry Smith Input Parameters: 4276b2b7091SBarry Smith + x - old solution vector 4286b2b7091SBarry Smith . y - search direction vector 4296b2b7091SBarry Smith . w - new solution vector 4306b2b7091SBarry Smith . changed_y - indicates that the line search changed y 4316b2b7091SBarry Smith - changed_w - indicates that the line search changed w 4326b2b7091SBarry Smith 433f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck() 434f190f2fcSBarry Smith and SNESLineSearchGetPostCheck() 435f190f2fcSBarry Smith 436878cb397SSatish Balay Level: advanced 4376b2b7091SBarry Smith 438f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck() 4396b2b7091SBarry Smith M*/ 4406b2b7091SBarry Smith 44186d74e61SPeter Brune #undef __FUNCT__ 442f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck" 44386d74e61SPeter Brune /*@C 444f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 445f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 44686d74e61SPeter Brune 447f1c6b773SPeter Brune Logically Collective on SNESLineSearch 44886d74e61SPeter Brune 44986d74e61SPeter Brune Input Parameters: 450f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 451f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence 452f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 45386d74e61SPeter Brune 45486d74e61SPeter Brune Level: intermediate 45586d74e61SPeter Brune 45686d74e61SPeter Brune .keywords: set, linesearch, post-check 45786d74e61SPeter Brune 458f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 45986d74e61SPeter Brune @*/ 460f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 46186d74e61SPeter Brune { 46286d74e61SPeter Brune PetscFunctionBegin; 463f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 464f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 46586d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 46686d74e61SPeter Brune PetscFunctionReturn(0); 46786d74e61SPeter Brune } 46886d74e61SPeter Brune 46986d74e61SPeter Brune #undef __FUNCT__ 470f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck" 47186d74e61SPeter Brune /*@C 472f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 47386d74e61SPeter Brune 47486d74e61SPeter Brune Input Parameters: 475f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 47686d74e61SPeter Brune 47786d74e61SPeter Brune Output Parameters: 478f190f2fcSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction 479f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 48086d74e61SPeter Brune 48186d74e61SPeter Brune Level: intermediate 48286d74e61SPeter Brune 48386d74e61SPeter Brune .keywords: get, linesearch, post-check 48486d74e61SPeter Brune 485f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck() 48686d74e61SPeter Brune @*/ 487f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 48886d74e61SPeter Brune { 4899bd66eb0SPeter Brune PetscFunctionBegin; 490f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 491f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 49286d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 49386d74e61SPeter Brune PetscFunctionReturn(0); 49486d74e61SPeter Brune } 49586d74e61SPeter Brune 496bf7f4e0aSPeter Brune #undef __FUNCT__ 497f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck" 498f40b411bSPeter Brune /*@ 499f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 500f40b411bSPeter Brune 501cd7522eaSPeter Brune Logically Collective on SNESLineSearch 502f40b411bSPeter Brune 503f40b411bSPeter Brune Input Parameters: 5047b1df9c1SPeter Brune + linesearch - The linesearch instance. 5057b1df9c1SPeter Brune . X - The current solution 5067b1df9c1SPeter Brune - Y - The step direction 507f40b411bSPeter Brune 508f40b411bSPeter Brune Output Parameters: 5098e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 510f40b411bSPeter Brune 511f190f2fcSBarry Smith Level: developer 512f40b411bSPeter Brune 513f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 514f40b411bSPeter Brune 515f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 516f40b411bSPeter Brune @*/ 5177b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 518bf7f4e0aSPeter Brune { 519bf7f4e0aSPeter Brune PetscErrorCode ierr; 5205fd66863SKarl Rupp 521bf7f4e0aSPeter Brune PetscFunctionBegin; 522bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 5236b2b7091SBarry Smith if (linesearch->ops->precheck) { 5246b2b7091SBarry Smith ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr); 52538bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 526bf7f4e0aSPeter Brune } 527bf7f4e0aSPeter Brune PetscFunctionReturn(0); 528bf7f4e0aSPeter Brune } 529bf7f4e0aSPeter Brune 530bf7f4e0aSPeter Brune #undef __FUNCT__ 531f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck" 532f40b411bSPeter Brune /*@ 533f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 534f40b411bSPeter Brune 535cd7522eaSPeter Brune Logically Collective on SNESLineSearch 536f40b411bSPeter Brune 537f40b411bSPeter Brune Input Parameters: 5387b1df9c1SPeter Brune + linesearch - The linesearch context 5397b1df9c1SPeter Brune . X - The last solution 5407b1df9c1SPeter Brune . Y - The step direction 5417b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 542f40b411bSPeter Brune 543f40b411bSPeter Brune Output Parameters: 54478bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 54578bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 546f40b411bSPeter Brune 547f190f2fcSBarry Smith Level: developer 548f40b411bSPeter Brune 549f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 550f40b411bSPeter Brune 551f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 552f40b411bSPeter Brune @*/ 5537b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 554bf7f4e0aSPeter Brune { 555bf7f4e0aSPeter Brune PetscErrorCode ierr; 556bf388a1fSBarry Smith 557bf7f4e0aSPeter Brune PetscFunctionBegin; 558bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 559bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 5606b2b7091SBarry Smith if (linesearch->ops->postcheck) { 5616b2b7091SBarry Smith ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr); 56238bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 56338bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 56486d74e61SPeter Brune } 56586d74e61SPeter Brune PetscFunctionReturn(0); 56686d74e61SPeter Brune } 56786d74e61SPeter Brune 56886d74e61SPeter Brune #undef __FUNCT__ 569f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard" 57086d74e61SPeter Brune /*@C 57186d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 57286d74e61SPeter Brune 573cd7522eaSPeter Brune Logically Collective on SNESLineSearch 57486d74e61SPeter Brune 57586d74e61SPeter Brune Input Arguments: 57686d74e61SPeter Brune + linesearch - linesearch context 57786d74e61SPeter Brune . X - base state for this step 57886d74e61SPeter Brune . Y - initial correction 579907376e6SBarry Smith - ctx - context for this function 58086d74e61SPeter Brune 58186d74e61SPeter Brune Output Arguments: 58286d74e61SPeter Brune + Y - correction, possibly modified 58386d74e61SPeter Brune - changed - flag indicating that Y was modified 58486d74e61SPeter Brune 58586d74e61SPeter Brune Options Database Key: 586cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 587cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 58886d74e61SPeter Brune 58986d74e61SPeter Brune Level: advanced 59086d74e61SPeter Brune 59186d74e61SPeter Brune Notes: 59286d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 59386d74e61SPeter Brune 59486d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 59586d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 59686d74e61SPeter Brune is generally not useful when using a Newton linearization. 59786d74e61SPeter Brune 59886d74e61SPeter Brune Reference: 59986d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 60086d74e61SPeter Brune 60186d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 60286d74e61SPeter Brune @*/ 603f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 60486d74e61SPeter Brune { 60586d74e61SPeter Brune PetscErrorCode ierr; 60686d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 60786d74e61SPeter Brune Vec Ylast; 60886d74e61SPeter Brune PetscScalar dot; 60986d74e61SPeter Brune PetscInt iter; 61086d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 61186d74e61SPeter Brune SNES snes; 61286d74e61SPeter Brune 61386d74e61SPeter Brune PetscFunctionBegin; 614f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 61586d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 61686d74e61SPeter Brune if (!Ylast) { 61786d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 61886d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 61986d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 62086d74e61SPeter Brune } 62186d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 62286d74e61SPeter Brune if (iter < 2) { 62386d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 62486d74e61SPeter Brune *changed = PETSC_FALSE; 62586d74e61SPeter Brune PetscFunctionReturn(0); 62686d74e61SPeter Brune } 62786d74e61SPeter Brune 62886d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 62986d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 63086d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 63186d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 632255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 63386d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 63486d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 63586d74e61SPeter Brune /* Modify the step Y */ 63686d74e61SPeter Brune PetscReal alpha,ydiffnorm; 63786d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 63886d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 63986d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 64086d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 64186d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 642c69d1a72SBarry 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); 64386d74e61SPeter Brune } else { 644c69d1a72SBarry 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); 64586d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 64686d74e61SPeter Brune *changed = PETSC_FALSE; 647bf7f4e0aSPeter Brune } 648bf7f4e0aSPeter Brune PetscFunctionReturn(0); 649bf7f4e0aSPeter Brune } 650bf7f4e0aSPeter Brune 651bf7f4e0aSPeter Brune #undef __FUNCT__ 652f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply" 653f40b411bSPeter Brune /*@ 654cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 655f40b411bSPeter Brune 656f1c6b773SPeter Brune Collective on SNESLineSearch 657f40b411bSPeter Brune 658f40b411bSPeter Brune Input Parameters: 6598e557f58SPeter Brune + linesearch - The linesearch context 6608e557f58SPeter Brune . X - The current solution 6618e557f58SPeter Brune . F - The current function 6628e557f58SPeter Brune . fnorm - The current norm 6638e557f58SPeter Brune - Y - The search direction 664f40b411bSPeter Brune 665f40b411bSPeter Brune Output Parameters: 6668e557f58SPeter Brune + X - The new solution 6678e557f58SPeter Brune . F - The new function 6688e557f58SPeter Brune - fnorm - The new function norm 669f40b411bSPeter Brune 670cd7522eaSPeter Brune Options Database Keys: 671d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell 672dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 6733c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 674cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 6753c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 6763c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 677cd7522eaSPeter Brune 678cd7522eaSPeter Brune Notes: 679cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 680cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 681cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 682cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 683cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 684cd7522eaSPeter Brune therefore may be fairly expensive. 685cd7522eaSPeter Brune 686f40b411bSPeter Brune Level: Intermediate 687f40b411bSPeter Brune 688f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 689f40b411bSPeter Brune 690cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction() 691f40b411bSPeter Brune @*/ 692bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 693bf388a1fSBarry Smith { 694bf7f4e0aSPeter Brune PetscErrorCode ierr; 695bf7f4e0aSPeter Brune 696bf388a1fSBarry Smith PetscFunctionBegin; 697f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 698bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 699bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 700bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 701bf7f4e0aSPeter Brune 702422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 703bf7f4e0aSPeter Brune 704bf7f4e0aSPeter Brune linesearch->vec_sol = X; 705bf7f4e0aSPeter Brune linesearch->vec_update = Y; 706bf7f4e0aSPeter Brune linesearch->vec_func = F; 707bf7f4e0aSPeter Brune 708f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 709bf7f4e0aSPeter Brune 710f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 711bf7f4e0aSPeter Brune 712f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 713f5af7f23SKarl Rupp else { 714bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 715bf7f4e0aSPeter Brune } 716bf7f4e0aSPeter Brune 717f1c6b773SPeter Brune ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 718bf7f4e0aSPeter Brune 719bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 720bf7f4e0aSPeter Brune 721f1c6b773SPeter Brune ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 722bf7f4e0aSPeter Brune 723f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 724bf7f4e0aSPeter Brune PetscFunctionReturn(0); 725bf7f4e0aSPeter Brune } 726bf7f4e0aSPeter Brune 727bf7f4e0aSPeter Brune #undef __FUNCT__ 728f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy" 729f40b411bSPeter Brune /*@ 730f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 731f40b411bSPeter Brune 732f1c6b773SPeter Brune Collective on SNESLineSearch 733f40b411bSPeter Brune 734f40b411bSPeter Brune Input Parameters: 7358e557f58SPeter Brune . linesearch - The linesearch context 736f40b411bSPeter Brune 737f40b411bSPeter Brune Level: Intermediate 738f40b411bSPeter Brune 73978bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 740f40b411bSPeter Brune 741cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 742f40b411bSPeter Brune @*/ 743bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 744bf388a1fSBarry Smith { 745bf7f4e0aSPeter Brune PetscErrorCode ierr; 746bf388a1fSBarry Smith 747bf7f4e0aSPeter Brune PetscFunctionBegin; 748bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 749f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 750bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 751e04113cfSBarry Smith ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr); 75222d28d08SBarry Smith ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr); 753f5af7f23SKarl Rupp if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch); 754bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 755dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr); 756e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 757bf7f4e0aSPeter Brune PetscFunctionReturn(0); 758bf7f4e0aSPeter Brune } 759bf7f4e0aSPeter Brune 760bf7f4e0aSPeter Brune #undef __FUNCT__ 761dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchSetDefaultMonitor" 762f40b411bSPeter Brune /*@ 763dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 764bf7f4e0aSPeter Brune 765bf7f4e0aSPeter Brune Input Parameters: 766dcf2fd19SBarry Smith + linesearch - the linesearch object 767dcf2fd19SBarry Smith - viewer - an ASCII PetscViewer or NULL to turn off monitor 768bf7f4e0aSPeter Brune 769dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 770bf7f4e0aSPeter Brune 771bf7f4e0aSPeter Brune Options Database: 772dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 773bf7f4e0aSPeter Brune 774bf7f4e0aSPeter Brune Level: intermediate 775bf7f4e0aSPeter Brune 776dcf2fd19SBarry Smith Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with 777d12e167eSBarry Smith SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the 778d12e167eSBarry Smith line search that are not visible to the other monitors. 779dcf2fd19SBarry Smith 780d12e167eSBarry Smith .seealso: SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor() 781bf7f4e0aSPeter Brune @*/ 782dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 783bf7f4e0aSPeter Brune { 784bf7f4e0aSPeter Brune PetscErrorCode ierr; 785bf388a1fSBarry Smith 786bf7f4e0aSPeter Brune PetscFunctionBegin; 787dcf2fd19SBarry Smith if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);} 788bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 789dcf2fd19SBarry Smith linesearch->monitor = viewer; 790bf7f4e0aSPeter Brune PetscFunctionReturn(0); 791bf7f4e0aSPeter Brune } 792bf7f4e0aSPeter Brune 793bf7f4e0aSPeter Brune #undef __FUNCT__ 794dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchGetDefaultMonitor" 795f40b411bSPeter Brune /*@ 796dcf2fd19SBarry Smith SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor. 7976a388c36SPeter Brune 798f190f2fcSBarry Smith Input Parameter: 7998e557f58SPeter Brune . linesearch - linesearch context 800f40b411bSPeter Brune 801f190f2fcSBarry Smith Output Parameter: 8028e557f58SPeter Brune . monitor - monitor context 803f40b411bSPeter Brune 804f40b411bSPeter Brune Logically Collective on SNES 805f40b411bSPeter Brune 806f40b411bSPeter Brune Options Database Keys: 8078e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 808f40b411bSPeter Brune 809f40b411bSPeter Brune Level: intermediate 810f40b411bSPeter Brune 811dcf2fd19SBarry Smith .seealso: SNESLineSearchSetDefaultMonitor(), PetscViewer 812f40b411bSPeter Brune @*/ 813dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 8146a388c36SPeter Brune { 8156a388c36SPeter Brune PetscFunctionBegin; 816f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8176a388c36SPeter Brune if (monitor) { 8186a388c36SPeter Brune PetscValidPointer(monitor, 2); 8196a388c36SPeter Brune *monitor = linesearch->monitor; 8206a388c36SPeter Brune } 8216a388c36SPeter Brune PetscFunctionReturn(0); 8226a388c36SPeter Brune } 8236a388c36SPeter Brune 8246a388c36SPeter Brune #undef __FUNCT__ 825dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSetFromOptions" 826dcf2fd19SBarry Smith /*@C 827dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 828dcf2fd19SBarry Smith 829dcf2fd19SBarry Smith Collective on SNESLineSearch 830dcf2fd19SBarry Smith 831dcf2fd19SBarry Smith Input Parameters: 832dcf2fd19SBarry Smith + ls - LineSearch object you wish to monitor 833dcf2fd19SBarry Smith . name - the monitor type one is seeking 834dcf2fd19SBarry Smith . help - message indicating what monitoring is done 835dcf2fd19SBarry Smith . manual - manual page for the monitor 836dcf2fd19SBarry Smith . monitor - the monitor function 837dcf2fd19SBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNESLineSearch or PetscViewer objects 838dcf2fd19SBarry Smith 839dcf2fd19SBarry Smith Level: developer 840dcf2fd19SBarry Smith 841dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 842dcf2fd19SBarry Smith PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 843dcf2fd19SBarry Smith PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 844dcf2fd19SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 845dcf2fd19SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 846dcf2fd19SBarry Smith PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 847dcf2fd19SBarry Smith PetscOptionsFList(), PetscOptionsEList() 848dcf2fd19SBarry Smith @*/ 849d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*)) 850dcf2fd19SBarry Smith { 851dcf2fd19SBarry Smith PetscErrorCode ierr; 852dcf2fd19SBarry Smith PetscViewer viewer; 853dcf2fd19SBarry Smith PetscViewerFormat format; 854dcf2fd19SBarry Smith PetscBool flg; 855dcf2fd19SBarry Smith 856dcf2fd19SBarry Smith PetscFunctionBegin; 857dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 858dcf2fd19SBarry Smith if (flg) { 859d12e167eSBarry Smith PetscViewerAndFormat *vf; 860d12e167eSBarry Smith ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 861d12e167eSBarry Smith ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 862dcf2fd19SBarry Smith if (monitorsetup) { 863d12e167eSBarry Smith ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr); 864dcf2fd19SBarry Smith } 865d12e167eSBarry Smith ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 866dcf2fd19SBarry Smith } 867dcf2fd19SBarry Smith PetscFunctionReturn(0); 868dcf2fd19SBarry Smith } 869dcf2fd19SBarry Smith 870dcf2fd19SBarry Smith #undef __FUNCT__ 871f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions" 872f40b411bSPeter Brune /*@ 873f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 874f40b411bSPeter Brune 875f40b411bSPeter Brune Input Parameters: 8768e557f58SPeter Brune . linesearch - linesearch context 877f40b411bSPeter Brune 878f40b411bSPeter Brune Options Database Keys: 879d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell 8803c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 8815a9b6599SPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch type 88271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 8831a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 8841a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 8851a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 8861a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 8871a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 888dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 889dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 8908e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 891cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 8921a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 8931a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method 894f40b411bSPeter Brune 895f1c6b773SPeter Brune Logically Collective on SNESLineSearch 896f40b411bSPeter Brune 897f40b411bSPeter Brune Level: intermediate 898f40b411bSPeter Brune 8993c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard() 900f40b411bSPeter Brune @*/ 901bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 902bf388a1fSBarry Smith { 903bf7f4e0aSPeter Brune PetscErrorCode ierr; 9041a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 905bf7f4e0aSPeter Brune char type[256]; 906bf7f4e0aSPeter Brune PetscBool flg, set; 907dcf2fd19SBarry Smith PetscViewer viewer; 908bf388a1fSBarry Smith 909bf7f4e0aSPeter Brune PetscFunctionBegin; 9100f51fdf8SToby Isaac ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr); 911bf7f4e0aSPeter Brune 912bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 913f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 914a264d7a6SBarry Smith ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 915bf7f4e0aSPeter Brune if (flg) { 916f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 917bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 918f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 919bf7f4e0aSPeter Brune } 920bf7f4e0aSPeter Brune 921dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr); 922dcf2fd19SBarry Smith if (set) { 923dcf2fd19SBarry Smith ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr); 924dcf2fd19SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 925dcf2fd19SBarry Smith } 926dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr); 927bf7f4e0aSPeter Brune 9281a9b3a06SPeter Brune /* tolerances */ 92994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr); 93094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr); 93194ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr); 93294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr); 93394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr); 93494ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr); 9357a35526eSPeter Brune 9361a9b3a06SPeter Brune /* damping parameters */ 93794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr); 9381a9b3a06SPeter Brune 93994ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr); 9401a9b3a06SPeter Brune 9411a9b3a06SPeter Brune /* precheck */ 9427a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 9437a35526eSPeter Brune if (set) { 9447a35526eSPeter Brune if (flg) { 9457a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 946f5af7f23SKarl Rupp 9477a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 9480298fd71SBarry Smith "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr); 9497a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 9507a35526eSPeter Brune } else { 9510298fd71SBarry Smith ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr); 9527a35526eSPeter Brune } 9537a35526eSPeter Brune } 95494ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr); 95594ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr); 9567a35526eSPeter Brune 9575a9b6599SPeter Brune if (linesearch->ops->setfromoptions) { 958e55864a3SBarry Smith (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr); 9595a9b6599SPeter Brune } 9605a9b6599SPeter Brune 9610633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr); 962bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 963bf7f4e0aSPeter Brune PetscFunctionReturn(0); 964bf7f4e0aSPeter Brune } 965bf7f4e0aSPeter Brune 966bf7f4e0aSPeter Brune #undef __FUNCT__ 967f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView" 968f40b411bSPeter Brune /*@ 969f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 970f40b411bSPeter Brune 971f40b411bSPeter Brune Input Parameters: 9728e557f58SPeter Brune . linesearch - linesearch context 973f40b411bSPeter Brune 974f1c6b773SPeter Brune Logically Collective on SNESLineSearch 975f40b411bSPeter Brune 976f40b411bSPeter Brune Level: intermediate 977f40b411bSPeter Brune 978dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate() 979f40b411bSPeter Brune @*/ 980bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 981bf388a1fSBarry Smith { 9827f1410a3SPeter Brune PetscErrorCode ierr; 9837f1410a3SPeter Brune PetscBool iascii; 984bf388a1fSBarry Smith 985bf7f4e0aSPeter Brune PetscFunctionBegin; 9867f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9877f1410a3SPeter Brune if (!viewer) { 988ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr); 9897f1410a3SPeter Brune } 9907f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 9917f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 992f40b411bSPeter Brune 993251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9947f1410a3SPeter Brune if (iascii) { 995dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr); 9967f1410a3SPeter Brune if (linesearch->ops->view) { 9977f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9987f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 9997f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 10007f1410a3SPeter Brune } 1001c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr); 1002c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr); 10037f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 10046b2b7091SBarry Smith if (linesearch->ops->precheck) { 10056b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 10067f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 10077f1410a3SPeter Brune } else { 10087f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 10097f1410a3SPeter Brune } 10107f1410a3SPeter Brune } 10116b2b7091SBarry Smith if (linesearch->ops->postcheck) { 10127f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 10137f1410a3SPeter Brune } 10147f1410a3SPeter Brune } 1015bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1016bf7f4e0aSPeter Brune } 1017bf7f4e0aSPeter Brune 1018bf7f4e0aSPeter Brune #undef __FUNCT__ 1019f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType" 1020ea5d4fccSPeter Brune /*@C 1021f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 1022f40b411bSPeter Brune 1023f190f2fcSBarry Smith Logically Collective on SNESLineSearch 1024f190f2fcSBarry Smith 1025f40b411bSPeter Brune Input Parameters: 10268e557f58SPeter Brune + linesearch - linesearch context 1027f40b411bSPeter Brune - type - The type of line search to be used 1028f40b411bSPeter Brune 1029cd7522eaSPeter Brune Available Types: 1030cd7522eaSPeter Brune + basic - Simple damping line search. 10318e557f58SPeter Brune . bt - Backtracking line search over the L2 norm of the function 10328e557f58SPeter Brune . l2 - Secant line search over the L2 norm of the function 10338e557f58SPeter Brune . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 1034d4c6564cSPatrick Farrell . nleqerr - Affine-covariant error-oriented linesearch 10358e557f58SPeter Brune - shell - User provided SNESLineSearch implementation 1036cd7522eaSPeter Brune 1037f40b411bSPeter Brune Level: intermediate 1038f40b411bSPeter Brune 1039f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 1040f40b411bSPeter Brune @*/ 104119fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 1042bf7f4e0aSPeter Brune { 1043f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 1044bf7f4e0aSPeter Brune PetscBool match; 1045bf7f4e0aSPeter Brune 1046bf7f4e0aSPeter Brune PetscFunctionBegin; 1047f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1048bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 1049bf7f4e0aSPeter Brune 1050251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 1051bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 1052bf7f4e0aSPeter Brune 10531c9cd337SJed Brown ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr); 1054bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 1055bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 1056bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 1057bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 1058f5af7f23SKarl Rupp 10590298fd71SBarry Smith linesearch->ops->destroy = NULL; 1060bf7f4e0aSPeter Brune } 1061f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 1062bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 1063bf7f4e0aSPeter Brune linesearch->ops->view = 0; 1064bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 1065bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 1066bf7f4e0aSPeter Brune 1067bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 1068bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 1069bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1070bf7f4e0aSPeter Brune } 1071bf7f4e0aSPeter Brune 1072bf7f4e0aSPeter Brune #undef __FUNCT__ 1073f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES" 1074f40b411bSPeter Brune /*@ 107578bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 1076f40b411bSPeter Brune 1077f40b411bSPeter Brune Input Parameters: 10788e557f58SPeter Brune + linesearch - linesearch context 1079f40b411bSPeter Brune - snes - The snes instance 1080f40b411bSPeter Brune 108178bcb3b5SPeter Brune Level: developer 108278bcb3b5SPeter Brune 108378bcb3b5SPeter Brune Notes: 1084f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 10857601faf0SJed Brown SNESGetLineSearch(). This routine is therefore mainly called within SNES 108678bcb3b5SPeter Brune implementations. 1087f40b411bSPeter Brune 10888141a3b9SPeter Brune Level: developer 1089f40b411bSPeter Brune 1090cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1091f40b411bSPeter Brune @*/ 1092bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1093bf388a1fSBarry Smith { 1094bf7f4e0aSPeter Brune PetscFunctionBegin; 1095f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1096bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 1097bf7f4e0aSPeter Brune linesearch->snes = snes; 1098bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1099bf7f4e0aSPeter Brune } 1100bf7f4e0aSPeter Brune 1101bf7f4e0aSPeter Brune #undef __FUNCT__ 1102f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES" 1103f40b411bSPeter Brune /*@ 11048141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 11058141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 11068141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 11078141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 1108f40b411bSPeter Brune 1109f40b411bSPeter Brune Input Parameters: 11108e557f58SPeter Brune . linesearch - linesearch context 1111f40b411bSPeter Brune 1112f40b411bSPeter Brune Output Parameters: 1113f40b411bSPeter Brune . snes - The snes instance 1114f40b411bSPeter Brune 11158141a3b9SPeter Brune Level: developer 1116f40b411bSPeter Brune 1117cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1118f40b411bSPeter Brune @*/ 1119bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1120bf388a1fSBarry Smith { 1121bf7f4e0aSPeter Brune PetscFunctionBegin; 1122f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11236a388c36SPeter Brune PetscValidPointer(snes, 2); 1124bf7f4e0aSPeter Brune *snes = linesearch->snes; 1125bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1126bf7f4e0aSPeter Brune } 1127bf7f4e0aSPeter Brune 11286a388c36SPeter Brune #undef __FUNCT__ 1129f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda" 1130f40b411bSPeter Brune /*@ 1131f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1132f40b411bSPeter Brune 1133f40b411bSPeter Brune Input Parameters: 11348e557f58SPeter Brune . linesearch - linesearch context 1135f40b411bSPeter Brune 1136f40b411bSPeter Brune Output Parameters: 1137cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 1138f40b411bSPeter Brune 113978bcb3b5SPeter Brune Level: advanced 114078bcb3b5SPeter Brune 11418e557f58SPeter Brune Notes: 11428e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 114378bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 114478bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 114578bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 114678bcb3b5SPeter Brune 1147cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 1148f40b411bSPeter Brune @*/ 1149f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 11506a388c36SPeter Brune { 11516a388c36SPeter Brune PetscFunctionBegin; 1152f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11536a388c36SPeter Brune PetscValidPointer(lambda, 2); 11546a388c36SPeter Brune *lambda = linesearch->lambda; 11556a388c36SPeter Brune PetscFunctionReturn(0); 11566a388c36SPeter Brune } 11576a388c36SPeter Brune 11586a388c36SPeter Brune #undef __FUNCT__ 1159f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda" 1160f40b411bSPeter Brune /*@ 1161f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 1162f40b411bSPeter Brune 1163f40b411bSPeter Brune Input Parameters: 11648e557f58SPeter Brune + linesearch - linesearch context 1165f40b411bSPeter Brune - lambda - The last steplength. 1166f40b411bSPeter Brune 1167cd7522eaSPeter Brune Notes: 1168f190f2fcSBarry Smith This routine is typically used within implementations of SNESLineSearchApply() 1169cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 1170cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1171cd7522eaSPeter Brune as an inner scaling parameter. 1172cd7522eaSPeter Brune 117378bcb3b5SPeter Brune Level: advanced 1174f40b411bSPeter Brune 1175f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 1176f40b411bSPeter Brune @*/ 1177f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 11786a388c36SPeter Brune { 11796a388c36SPeter Brune PetscFunctionBegin; 1180f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11816a388c36SPeter Brune linesearch->lambda = lambda; 11826a388c36SPeter Brune PetscFunctionReturn(0); 11836a388c36SPeter Brune } 11846a388c36SPeter Brune 11856a388c36SPeter Brune #undef __FUNCT__ 1186f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances" 1187f40b411bSPeter Brune /*@ 11883c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 118978bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 119078bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 119178bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1192f40b411bSPeter Brune 1193f40b411bSPeter Brune Input Parameters: 11948e557f58SPeter Brune . linesearch - linesearch context 1195f40b411bSPeter Brune 1196f40b411bSPeter Brune Output Parameters: 1197516fe3c3SPeter Brune + steptol - The minimum steplength 11986cc8e53bSPeter Brune . maxstep - The maximum steplength 1199516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1200516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1201516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1202516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1203f40b411bSPeter Brune 120478bcb3b5SPeter Brune Level: intermediate 120578bcb3b5SPeter Brune 120678bcb3b5SPeter Brune Notes: 120778bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 12083c7d6663SPeter Brune the type requires. 1209516fe3c3SPeter Brune 1210f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 1211f40b411bSPeter Brune @*/ 1212f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 12136a388c36SPeter Brune { 12146a388c36SPeter Brune PetscFunctionBegin; 1215f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1216516fe3c3SPeter Brune if (steptol) { 12176a388c36SPeter Brune PetscValidPointer(steptol, 2); 12186a388c36SPeter Brune *steptol = linesearch->steptol; 1219516fe3c3SPeter Brune } 1220516fe3c3SPeter Brune if (maxstep) { 1221516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 1222516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1223516fe3c3SPeter Brune } 1224516fe3c3SPeter Brune if (rtol) { 1225516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 1226516fe3c3SPeter Brune *rtol = linesearch->rtol; 1227516fe3c3SPeter Brune } 1228516fe3c3SPeter Brune if (atol) { 1229516fe3c3SPeter Brune PetscValidPointer(atol, 5); 1230516fe3c3SPeter Brune *atol = linesearch->atol; 1231516fe3c3SPeter Brune } 1232516fe3c3SPeter Brune if (ltol) { 1233516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 1234516fe3c3SPeter Brune *ltol = linesearch->ltol; 1235516fe3c3SPeter Brune } 1236516fe3c3SPeter Brune if (max_its) { 1237516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 1238516fe3c3SPeter Brune *max_its = linesearch->max_its; 1239516fe3c3SPeter Brune } 12406a388c36SPeter Brune PetscFunctionReturn(0); 12416a388c36SPeter Brune } 12426a388c36SPeter Brune 12436a388c36SPeter Brune #undef __FUNCT__ 1244f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances" 1245f40b411bSPeter Brune /*@ 12463c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 124778bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 124878bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 124978bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1250f40b411bSPeter Brune 1251f40b411bSPeter Brune Input Parameters: 12528e557f58SPeter Brune + linesearch - linesearch context 1253516fe3c3SPeter Brune . steptol - The minimum steplength 12546cc8e53bSPeter Brune . maxstep - The maximum steplength 1255516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1256516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1257516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1258516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1259f40b411bSPeter Brune 126078bcb3b5SPeter Brune Notes: 12613c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 126278bcb3b5SPeter Brune place of an argument. 1263f40b411bSPeter Brune 126478bcb3b5SPeter Brune Level: intermediate 1265516fe3c3SPeter Brune 1266f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 1267f40b411bSPeter Brune @*/ 1268f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 12696a388c36SPeter Brune { 12706a388c36SPeter Brune PetscFunctionBegin; 1271f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1272d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1273d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1274d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1275d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1276d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1277d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1278d3952184SSatish Balay 1279d3952184SSatish Balay if (steptol!= PETSC_DEFAULT) { 1280ce94432eSBarry Smith if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 12816a388c36SPeter Brune linesearch->steptol = steptol; 1282d3952184SSatish Balay } 1283d3952184SSatish Balay 1284d3952184SSatish Balay if (maxstep!= PETSC_DEFAULT) { 1285ce94432eSBarry Smith if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1286516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1287d3952184SSatish Balay } 1288d3952184SSatish Balay 1289d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1290ce94432eSBarry Smith if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol); 1291516fe3c3SPeter Brune linesearch->rtol = rtol; 1292d3952184SSatish Balay } 1293d3952184SSatish Balay 1294d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1295ce94432eSBarry Smith if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1296516fe3c3SPeter Brune linesearch->atol = atol; 1297d3952184SSatish Balay } 1298d3952184SSatish Balay 1299d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1300ce94432eSBarry Smith if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1301516fe3c3SPeter Brune linesearch->ltol = ltol; 1302d3952184SSatish Balay } 1303d3952184SSatish Balay 1304d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1305ce94432eSBarry Smith if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1306516fe3c3SPeter Brune linesearch->max_its = max_its; 1307d3952184SSatish Balay } 13086a388c36SPeter Brune PetscFunctionReturn(0); 13096a388c36SPeter Brune } 13106a388c36SPeter Brune 13116a388c36SPeter Brune #undef __FUNCT__ 1312f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping" 1313f40b411bSPeter Brune /*@ 1314f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1315f40b411bSPeter Brune 1316f40b411bSPeter Brune Input Parameters: 13178e557f58SPeter Brune . linesearch - linesearch context 1318f40b411bSPeter Brune 1319f40b411bSPeter Brune Output Parameters: 13208e557f58SPeter Brune . damping - The damping parameter 1321f40b411bSPeter Brune 132278bcb3b5SPeter Brune Level: advanced 1323f40b411bSPeter Brune 132478bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1325f40b411bSPeter Brune @*/ 1326f40b411bSPeter Brune 1327f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 13286a388c36SPeter Brune { 13296a388c36SPeter Brune PetscFunctionBegin; 1330f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13316a388c36SPeter Brune PetscValidPointer(damping, 2); 13326a388c36SPeter Brune *damping = linesearch->damping; 13336a388c36SPeter Brune PetscFunctionReturn(0); 13346a388c36SPeter Brune } 13356a388c36SPeter Brune 13366a388c36SPeter Brune #undef __FUNCT__ 1337f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping" 1338f40b411bSPeter Brune /*@ 1339f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1340f40b411bSPeter Brune 1341f40b411bSPeter Brune Input Parameters: 134203fd4120SBarry Smith + linesearch - linesearch context 134303fd4120SBarry Smith - damping - The damping parameter 1344f40b411bSPeter Brune 134503fd4120SBarry Smith Options Database: 134603fd4120SBarry Smith . -snes_linesearch_damping 1347f40b411bSPeter Brune Level: intermediate 1348f40b411bSPeter Brune 1349cd7522eaSPeter Brune Notes: 1350cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1351cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 135278bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1353cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1354cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1355cd7522eaSPeter Brune 1356f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1357f40b411bSPeter Brune @*/ 1358f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 13596a388c36SPeter Brune { 13606a388c36SPeter Brune PetscFunctionBegin; 1361f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13626a388c36SPeter Brune linesearch->damping = damping; 13636a388c36SPeter Brune PetscFunctionReturn(0); 13646a388c36SPeter Brune } 13656a388c36SPeter Brune 13666a388c36SPeter Brune #undef __FUNCT__ 136759405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder" 136859405d5eSPeter Brune /*@ 136959405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 137059405d5eSPeter Brune 137159405d5eSPeter Brune Input Parameters: 137278bcb3b5SPeter Brune . linesearch - linesearch context 137359405d5eSPeter Brune 137459405d5eSPeter Brune Output Parameters: 13758e557f58SPeter Brune . order - The order 137659405d5eSPeter Brune 137778bcb3b5SPeter Brune Possible Values for order: 13783c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 13793c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 13803c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 138178bcb3b5SPeter Brune 138259405d5eSPeter Brune Level: intermediate 138359405d5eSPeter Brune 138459405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 138559405d5eSPeter Brune @*/ 138659405d5eSPeter Brune 1387b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 138859405d5eSPeter Brune { 138959405d5eSPeter Brune PetscFunctionBegin; 139059405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 139159405d5eSPeter Brune PetscValidPointer(order, 2); 139259405d5eSPeter Brune *order = linesearch->order; 139359405d5eSPeter Brune PetscFunctionReturn(0); 139459405d5eSPeter Brune } 139559405d5eSPeter Brune 139659405d5eSPeter Brune #undef __FUNCT__ 139759405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder" 139859405d5eSPeter Brune /*@ 139959405d5eSPeter Brune SNESLineSearchSetOrder - Sets the line search damping paramter. 140059405d5eSPeter Brune 140159405d5eSPeter Brune Input Parameters: 140278bcb3b5SPeter Brune . linesearch - linesearch context 140378bcb3b5SPeter Brune . order - The damping parameter 140459405d5eSPeter Brune 140559405d5eSPeter Brune Level: intermediate 140659405d5eSPeter Brune 140778bcb3b5SPeter Brune Possible Values for order: 14083c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 14093c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 14103c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 141178bcb3b5SPeter Brune 141259405d5eSPeter Brune Notes: 141359405d5eSPeter Brune Variable orders are supported by the following line searches: 141478bcb3b5SPeter Brune + bt - cubic and quadratic 141578bcb3b5SPeter Brune - cp - linear and quadratic 141659405d5eSPeter Brune 141759405d5eSPeter Brune .seealso: SNESLineSearchGetOrder() 141859405d5eSPeter Brune @*/ 1419b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 142059405d5eSPeter Brune { 142159405d5eSPeter Brune PetscFunctionBegin; 142259405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 142359405d5eSPeter Brune linesearch->order = order; 142459405d5eSPeter Brune PetscFunctionReturn(0); 142559405d5eSPeter Brune } 142659405d5eSPeter Brune 142759405d5eSPeter Brune #undef __FUNCT__ 1428f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms" 1429f40b411bSPeter Brune /*@ 1430f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1431f40b411bSPeter Brune 1432f40b411bSPeter Brune Input Parameters: 143378bcb3b5SPeter Brune . linesearch - linesearch context 1434f40b411bSPeter Brune 1435f40b411bSPeter Brune Output Parameters: 1436f40b411bSPeter Brune + xnorm - The norm of the current solution 1437f40b411bSPeter Brune . fnorm - The norm of the current function 1438f40b411bSPeter Brune - ynorm - The norm of the current update 1439f40b411bSPeter Brune 1440cd7522eaSPeter Brune Notes: 1441cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1442cd7522eaSPeter Brune 144378bcb3b5SPeter Brune Level: developer 1444f40b411bSPeter Brune 1445f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1446f40b411bSPeter Brune @*/ 1447f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1448bf7f4e0aSPeter Brune { 1449bf7f4e0aSPeter Brune PetscFunctionBegin; 1450f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1451f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1452f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1453f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 1454bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1455bf7f4e0aSPeter Brune } 1456bf7f4e0aSPeter Brune 1457e7058c64SPeter Brune #undef __FUNCT__ 1458f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms" 1459f40b411bSPeter Brune /*@ 1460f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1461f40b411bSPeter Brune 1462f40b411bSPeter Brune Input Parameters: 146378bcb3b5SPeter Brune + linesearch - linesearch context 1464f40b411bSPeter Brune . xnorm - The norm of the current solution 1465f40b411bSPeter Brune . fnorm - The norm of the current function 1466f40b411bSPeter Brune - ynorm - The norm of the current update 1467f40b411bSPeter Brune 146878bcb3b5SPeter Brune Level: advanced 1469f40b411bSPeter Brune 1470f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1471f40b411bSPeter Brune @*/ 1472f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 14736a388c36SPeter Brune { 14746a388c36SPeter Brune PetscFunctionBegin; 1475f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14766a388c36SPeter Brune linesearch->xnorm = xnorm; 14776a388c36SPeter Brune linesearch->fnorm = fnorm; 14786a388c36SPeter Brune linesearch->ynorm = ynorm; 14796a388c36SPeter Brune PetscFunctionReturn(0); 14806a388c36SPeter Brune } 14816a388c36SPeter Brune 14826a388c36SPeter Brune #undef __FUNCT__ 1483f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms" 1484f40b411bSPeter Brune /*@ 1485f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1486f40b411bSPeter Brune 1487f40b411bSPeter Brune Input Parameters: 148878bcb3b5SPeter Brune . linesearch - linesearch context 1489f40b411bSPeter Brune 1490f40b411bSPeter Brune Options Database Keys: 14918e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1492f40b411bSPeter Brune 1493f40b411bSPeter Brune Level: intermediate 1494f40b411bSPeter Brune 149578bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1496f40b411bSPeter Brune @*/ 1497f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 14986a388c36SPeter Brune { 14996a388c36SPeter Brune PetscErrorCode ierr; 15009bd66eb0SPeter Brune SNES snes; 1501bf388a1fSBarry Smith 15026a388c36SPeter Brune PetscFunctionBegin; 15036a388c36SPeter Brune if (linesearch->norms) { 15049bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1505f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 15069bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 15079bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15089bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 15099bd66eb0SPeter Brune } else { 15106a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 15116a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 15126a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15136a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 15146a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 15156a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15166a388c36SPeter Brune } 15179bd66eb0SPeter Brune } 15186a388c36SPeter Brune PetscFunctionReturn(0); 15196a388c36SPeter Brune } 15206a388c36SPeter Brune 15216f263ca3SPeter Brune #undef __FUNCT__ 15226f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms" 15236f263ca3SPeter Brune /*@ 15246f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 15256f263ca3SPeter Brune 15266f263ca3SPeter Brune Input Parameters: 152778bcb3b5SPeter Brune + linesearch - linesearch context 152878bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 15296f263ca3SPeter Brune 15306f263ca3SPeter Brune Options Database Keys: 15318e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 15326f263ca3SPeter Brune 15336f263ca3SPeter Brune Notes: 15341a4f838cSPeter Brune This is most relevant to the SNESLINESEARCHBASIC line search type. 15356f263ca3SPeter Brune 15366f263ca3SPeter Brune Level: intermediate 15376f263ca3SPeter Brune 15381a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 15396f263ca3SPeter Brune @*/ 15406f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 15416f263ca3SPeter Brune { 15426f263ca3SPeter Brune PetscFunctionBegin; 15436f263ca3SPeter Brune linesearch->norms = flg; 15446f263ca3SPeter Brune PetscFunctionReturn(0); 15456f263ca3SPeter Brune } 15466f263ca3SPeter Brune 15476a388c36SPeter Brune #undef __FUNCT__ 1548f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs" 1549f40b411bSPeter Brune /*@ 1550f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1551f40b411bSPeter Brune 1552f40b411bSPeter Brune Input Parameters: 155378bcb3b5SPeter Brune . linesearch - linesearch context 1554f40b411bSPeter Brune 1555f40b411bSPeter Brune Output Parameters: 15566232e825SPeter Brune + X - Solution vector 15576232e825SPeter Brune . F - Function vector 15586232e825SPeter Brune . Y - Search direction vector 15596232e825SPeter Brune . W - Solution work vector 15606232e825SPeter Brune - G - Function work vector 15616232e825SPeter Brune 15627bba9028SPeter Brune Notes: 15637bba9028SPeter Brune At the beginning of a line search application, X should contain a 15646232e825SPeter Brune solution and the vector F the function computed at X. At the end of the 15656232e825SPeter Brune line search application, X should contain the new solution, and F the 15666232e825SPeter Brune function evaluated at the new solution. 1567f40b411bSPeter Brune 15682a7a6963SBarry Smith These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 15692a7a6963SBarry Smith 157078bcb3b5SPeter Brune Level: advanced 1571f40b411bSPeter Brune 1572f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1573f40b411bSPeter Brune @*/ 1574bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1575bf388a1fSBarry Smith { 15766a388c36SPeter Brune PetscFunctionBegin; 1577f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15786a388c36SPeter Brune if (X) { 15796a388c36SPeter Brune PetscValidPointer(X, 2); 15806a388c36SPeter Brune *X = linesearch->vec_sol; 15816a388c36SPeter Brune } 15826a388c36SPeter Brune if (F) { 15836a388c36SPeter Brune PetscValidPointer(F, 3); 15846a388c36SPeter Brune *F = linesearch->vec_func; 15856a388c36SPeter Brune } 15866a388c36SPeter Brune if (Y) { 15876a388c36SPeter Brune PetscValidPointer(Y, 4); 15886a388c36SPeter Brune *Y = linesearch->vec_update; 15896a388c36SPeter Brune } 15906a388c36SPeter Brune if (W) { 15916a388c36SPeter Brune PetscValidPointer(W, 5); 15926a388c36SPeter Brune *W = linesearch->vec_sol_new; 15936a388c36SPeter Brune } 15946a388c36SPeter Brune if (G) { 15956a388c36SPeter Brune PetscValidPointer(G, 6); 15966a388c36SPeter Brune *G = linesearch->vec_func_new; 15976a388c36SPeter Brune } 15986a388c36SPeter Brune PetscFunctionReturn(0); 15996a388c36SPeter Brune } 16006a388c36SPeter Brune 16016a388c36SPeter Brune #undef __FUNCT__ 1602f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs" 1603f40b411bSPeter Brune /*@ 1604f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1605f40b411bSPeter Brune 1606f40b411bSPeter Brune Input Parameters: 160778bcb3b5SPeter Brune + linesearch - linesearch context 16086232e825SPeter Brune . X - Solution vector 16096232e825SPeter Brune . F - Function vector 16106232e825SPeter Brune . Y - Search direction vector 16116232e825SPeter Brune . W - Solution work vector 16126232e825SPeter Brune - G - Function work vector 1613f40b411bSPeter Brune 161478bcb3b5SPeter Brune Level: advanced 1615f40b411bSPeter Brune 1616f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1617f40b411bSPeter Brune @*/ 1618bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1619bf388a1fSBarry Smith { 16206a388c36SPeter Brune PetscFunctionBegin; 1621f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 16226a388c36SPeter Brune if (X) { 16236a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 16246a388c36SPeter Brune linesearch->vec_sol = X; 16256a388c36SPeter Brune } 16266a388c36SPeter Brune if (F) { 16276a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 16286a388c36SPeter Brune linesearch->vec_func = F; 16296a388c36SPeter Brune } 16306a388c36SPeter Brune if (Y) { 16316a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 16326a388c36SPeter Brune linesearch->vec_update = Y; 16336a388c36SPeter Brune } 16346a388c36SPeter Brune if (W) { 16356a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 16366a388c36SPeter Brune linesearch->vec_sol_new = W; 16376a388c36SPeter Brune } 16386a388c36SPeter Brune if (G) { 16396a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 16406a388c36SPeter Brune linesearch->vec_func_new = G; 16416a388c36SPeter Brune } 16426a388c36SPeter Brune PetscFunctionReturn(0); 16436a388c36SPeter Brune } 16446a388c36SPeter Brune 16456a388c36SPeter Brune #undef __FUNCT__ 1646f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1647e7058c64SPeter Brune /*@C 1648f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1649e7058c64SPeter Brune SNES options in the database. 1650e7058c64SPeter Brune 1651cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1652e7058c64SPeter Brune 1653e7058c64SPeter Brune Input Parameters: 1654e7058c64SPeter Brune + snes - the SNES context 1655e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1656e7058c64SPeter Brune 1657e7058c64SPeter Brune Notes: 1658e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1659e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1660e7058c64SPeter Brune 1661e7058c64SPeter Brune Level: advanced 1662e7058c64SPeter Brune 1663f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1664e7058c64SPeter Brune 1665e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1666e7058c64SPeter Brune @*/ 1667f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1668e7058c64SPeter Brune { 1669e7058c64SPeter Brune PetscErrorCode ierr; 1670e7058c64SPeter Brune 1671e7058c64SPeter Brune PetscFunctionBegin; 1672f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1673e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1674e7058c64SPeter Brune PetscFunctionReturn(0); 1675e7058c64SPeter Brune } 1676e7058c64SPeter Brune 1677e7058c64SPeter Brune #undef __FUNCT__ 1678f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1679e7058c64SPeter Brune /*@C 1680f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1681f1c6b773SPeter Brune SNESLineSearch options in the database. 1682e7058c64SPeter Brune 1683e7058c64SPeter Brune Not Collective 1684e7058c64SPeter Brune 1685e7058c64SPeter Brune Input Parameter: 1686cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1687e7058c64SPeter Brune 1688e7058c64SPeter Brune Output Parameter: 1689e7058c64SPeter Brune . prefix - pointer to the prefix string used 1690e7058c64SPeter Brune 16918e557f58SPeter Brune Notes: 16928e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1693e7058c64SPeter Brune sufficient length to hold the prefix. 1694e7058c64SPeter Brune 1695e7058c64SPeter Brune Level: advanced 1696e7058c64SPeter Brune 1697f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1698e7058c64SPeter Brune 1699e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1700e7058c64SPeter Brune @*/ 1701f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1702e7058c64SPeter Brune { 1703e7058c64SPeter Brune PetscErrorCode ierr; 1704e7058c64SPeter Brune 1705e7058c64SPeter Brune PetscFunctionBegin; 1706f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1707e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1708e7058c64SPeter Brune PetscFunctionReturn(0); 1709e7058c64SPeter Brune } 1710bf7f4e0aSPeter Brune 1711bf7f4e0aSPeter Brune #undef __FUNCT__ 1712fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs" 17138d359177SBarry Smith /*@C 1714fa0ddf94SBarry Smith SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1715f40b411bSPeter Brune 1716f40b411bSPeter Brune Input Parameter: 1717f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1718f40b411bSPeter Brune - nwork - the number of work vectors 1719f40b411bSPeter Brune 1720f40b411bSPeter Brune Level: developer 1721f40b411bSPeter Brune 1722fa0ddf94SBarry Smith Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1723cd7522eaSPeter Brune 1724f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1725f40b411bSPeter Brune 1726fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs() 1727f40b411bSPeter Brune @*/ 1728fa0ddf94SBarry Smith PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1729bf7f4e0aSPeter Brune { 1730bf7f4e0aSPeter Brune PetscErrorCode ierr; 1731bf388a1fSBarry Smith 1732bf7f4e0aSPeter Brune PetscFunctionBegin; 1733bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1734bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 17358d359177SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1736bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1737bf7f4e0aSPeter Brune } 1738bf7f4e0aSPeter Brune 1739bf7f4e0aSPeter Brune #undef __FUNCT__ 1740422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchGetReason" 1741f40b411bSPeter Brune /*@ 1742422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1743f40b411bSPeter Brune 1744f40b411bSPeter Brune Input Parameters: 174578bcb3b5SPeter Brune . linesearch - linesearch context 1746f40b411bSPeter Brune 1747f40b411bSPeter Brune Output Parameters: 1748422a814eSBarry Smith . result - The success or failure status 1749f40b411bSPeter Brune 1750cd7522eaSPeter Brune Notes: 1751c98378a5SBarry Smith This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1752cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1753cd7522eaSPeter Brune 1754f40b411bSPeter Brune Level: intermediate 1755f40b411bSPeter Brune 1756422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason 1757f40b411bSPeter Brune @*/ 1758422a814eSBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1759bf7f4e0aSPeter Brune { 1760bf7f4e0aSPeter Brune PetscFunctionBegin; 1761f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1762422a814eSBarry Smith PetscValidPointer(result, 2); 1763422a814eSBarry Smith *result = linesearch->result; 1764bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1765bf7f4e0aSPeter Brune } 1766bf7f4e0aSPeter Brune 1767bf7f4e0aSPeter Brune #undef __FUNCT__ 1768422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchSetReason" 1769f40b411bSPeter Brune /*@ 1770422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1771f40b411bSPeter Brune 1772f40b411bSPeter Brune Input Parameters: 177378bcb3b5SPeter Brune + linesearch - linesearch context 1774422a814eSBarry Smith - result - The success or failure status 1775f40b411bSPeter Brune 1776cd7522eaSPeter Brune Notes: 1777cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1778cd7522eaSPeter Brune the success or failure of the line search method. 1779cd7522eaSPeter Brune 178078bcb3b5SPeter Brune Level: developer 1781f40b411bSPeter Brune 1782422a814eSBarry Smith .seealso: SNESLineSearchGetSResult() 1783f40b411bSPeter Brune @*/ 1784422a814eSBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 17856a388c36SPeter Brune { 17866a388c36SPeter Brune PetscFunctionBegin; 17875fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1788422a814eSBarry Smith linesearch->result = result; 17896a388c36SPeter Brune PetscFunctionReturn(0); 17906a388c36SPeter Brune } 17916a388c36SPeter Brune 17926a388c36SPeter Brune #undef __FUNCT__ 1793f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions" 17949bd66eb0SPeter Brune /*@C 1795f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 17969bd66eb0SPeter Brune 17979bd66eb0SPeter Brune Input Parameters: 17989bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 17999bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 18009bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 18019bd66eb0SPeter Brune 18029bd66eb0SPeter Brune Logically Collective on SNES 18039bd66eb0SPeter Brune 18049bd66eb0SPeter Brune Calling sequence of projectfunc: 18059bd66eb0SPeter Brune .vb 18069bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 18079bd66eb0SPeter Brune .ve 18089bd66eb0SPeter Brune 18099bd66eb0SPeter Brune Input parameters for projectfunc: 18109bd66eb0SPeter Brune + snes - nonlinear context 18119bd66eb0SPeter Brune - X - current solution 18129bd66eb0SPeter Brune 1813cd7522eaSPeter Brune Output parameters for projectfunc: 18149bd66eb0SPeter Brune . X - Projected solution 18159bd66eb0SPeter Brune 18169bd66eb0SPeter Brune Calling sequence of normfunc: 18179bd66eb0SPeter Brune .vb 18189bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 18199bd66eb0SPeter Brune .ve 18209bd66eb0SPeter Brune 1821cd7522eaSPeter Brune Input parameters for normfunc: 18229bd66eb0SPeter Brune + snes - nonlinear context 18239bd66eb0SPeter Brune . X - current solution 18249bd66eb0SPeter Brune - F - current residual 18259bd66eb0SPeter Brune 1826cd7522eaSPeter Brune Output parameters for normfunc: 18279bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 18289bd66eb0SPeter Brune 1829cd7522eaSPeter Brune Notes: 1830cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1831cd7522eaSPeter Brune 1832cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1833cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 18349bd66eb0SPeter Brune 18359bd66eb0SPeter Brune Level: developer 18369bd66eb0SPeter Brune 18379bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 18389bd66eb0SPeter Brune 1839f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 18409bd66eb0SPeter Brune @*/ 1841f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 18429bd66eb0SPeter Brune { 18439bd66eb0SPeter Brune PetscFunctionBegin; 1844f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 18459bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 18469bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 18479bd66eb0SPeter Brune PetscFunctionReturn(0); 18489bd66eb0SPeter Brune } 18499bd66eb0SPeter Brune 18509bd66eb0SPeter Brune #undef __FUNCT__ 1851f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions" 18529bd66eb0SPeter Brune /*@C 1853f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 18549bd66eb0SPeter Brune 18559bd66eb0SPeter Brune Input Parameters: 1856907376e6SBarry Smith . linesearch - the line search context, obtain with SNESGetLineSearch() 18579bd66eb0SPeter Brune 18589bd66eb0SPeter Brune Output Parameters: 18599bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 18609bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 18619bd66eb0SPeter Brune 18629bd66eb0SPeter Brune Logically Collective on SNES 18639bd66eb0SPeter Brune 18649bd66eb0SPeter Brune Level: developer 18659bd66eb0SPeter Brune 18669bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 18679bd66eb0SPeter Brune 1868f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 18699bd66eb0SPeter Brune @*/ 1870f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 18719bd66eb0SPeter Brune { 18729bd66eb0SPeter Brune PetscFunctionBegin; 18739bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 18749bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 18759bd66eb0SPeter Brune PetscFunctionReturn(0); 18769bd66eb0SPeter Brune } 18779bd66eb0SPeter Brune 18789bd66eb0SPeter Brune #undef __FUNCT__ 1879f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister" 1880bf7f4e0aSPeter Brune /*@C 18811c84c290SBarry Smith SNESLineSearchRegister - See SNESLineSearchRegister() 1882bf7f4e0aSPeter Brune 1883bf7f4e0aSPeter Brune Level: advanced 1884bf7f4e0aSPeter Brune @*/ 1885bdf89e91SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1886bf7f4e0aSPeter Brune { 1887bf7f4e0aSPeter Brune PetscErrorCode ierr; 1888bf7f4e0aSPeter Brune 1889bf7f4e0aSPeter Brune PetscFunctionBegin; 1890a240a19fSJed Brown ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr); 1891bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1892bf7f4e0aSPeter Brune } 1893