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; 757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply; 8bf7f4e0aSPeter Brune 9dcf2fd19SBarry Smith /*@ 10dcf2fd19SBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object. 11dcf2fd19SBarry Smith 12dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 13dcf2fd19SBarry Smith 14dcf2fd19SBarry Smith Input Parameters: 15dcf2fd19SBarry Smith . ls - the SNESLineSearch context 16dcf2fd19SBarry Smith 17dcf2fd19SBarry Smith Options Database Key: 18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 19dcf2fd19SBarry Smith into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those 20dcf2fd19SBarry Smith set via the options database 21dcf2fd19SBarry Smith 22dcf2fd19SBarry Smith Notes: 23dcf2fd19SBarry Smith There is no way to clear one specific monitor from a SNESLineSearch object. 24dcf2fd19SBarry Smith 25dcf2fd19SBarry Smith This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel 26dcf2fd19SBarry Smith that one. 27dcf2fd19SBarry Smith 28dcf2fd19SBarry Smith Level: intermediate 29dcf2fd19SBarry Smith 30dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 31dcf2fd19SBarry Smith 32dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet() 33dcf2fd19SBarry Smith @*/ 34dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 35dcf2fd19SBarry Smith { 36dcf2fd19SBarry Smith PetscErrorCode ierr; 37dcf2fd19SBarry Smith PetscInt i; 38dcf2fd19SBarry Smith 39dcf2fd19SBarry Smith PetscFunctionBegin; 40dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 41dcf2fd19SBarry Smith for (i=0; i<ls->numbermonitors; i++) { 42dcf2fd19SBarry Smith if (ls->monitordestroy[i]) { 43dcf2fd19SBarry Smith ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr); 44dcf2fd19SBarry Smith } 45dcf2fd19SBarry Smith } 46dcf2fd19SBarry Smith ls->numbermonitors = 0; 47dcf2fd19SBarry Smith PetscFunctionReturn(0); 48dcf2fd19SBarry Smith } 49dcf2fd19SBarry Smith 50dcf2fd19SBarry Smith /*@ 51dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 52dcf2fd19SBarry Smith 53dcf2fd19SBarry Smith Collective on SNES 54dcf2fd19SBarry Smith 55dcf2fd19SBarry Smith Input Parameters: 56dcf2fd19SBarry Smith . ls - the linesearch object 57dcf2fd19SBarry Smith 58dcf2fd19SBarry Smith Notes: 59dcf2fd19SBarry Smith This routine is called by the SNES implementations. 60dcf2fd19SBarry Smith It does not typically need to be called by the user. 61dcf2fd19SBarry Smith 62dcf2fd19SBarry Smith Level: developer 63dcf2fd19SBarry Smith 64dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorSet() 65dcf2fd19SBarry Smith @*/ 66dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 67dcf2fd19SBarry Smith { 68dcf2fd19SBarry Smith PetscErrorCode ierr; 69dcf2fd19SBarry Smith PetscInt i,n = ls->numbermonitors; 70dcf2fd19SBarry Smith 71dcf2fd19SBarry Smith PetscFunctionBegin; 72dcf2fd19SBarry Smith for (i=0; i<n; i++) { 73dcf2fd19SBarry Smith ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr); 74dcf2fd19SBarry Smith } 75dcf2fd19SBarry Smith PetscFunctionReturn(0); 76dcf2fd19SBarry Smith } 77dcf2fd19SBarry Smith 78dcf2fd19SBarry Smith /*@C 79dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 80dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 81dcf2fd19SBarry Smith progress. 82dcf2fd19SBarry Smith 83dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 84dcf2fd19SBarry Smith 85dcf2fd19SBarry Smith Input Parameters: 86dcf2fd19SBarry Smith + ls - the SNESLineSearch context 87dcf2fd19SBarry Smith . f - the monitor function 88dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 89dcf2fd19SBarry Smith monitor routine (use NULL if no context is desired) 90dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 91dcf2fd19SBarry Smith (may be NULL) 92dcf2fd19SBarry Smith 93dcf2fd19SBarry Smith Notes: 94dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 95dcf2fd19SBarry Smith SNESLineSearchMonitorSet() multiple times; all will be called in the 96dcf2fd19SBarry Smith order in which they were set. 97dcf2fd19SBarry Smith 98*95452b02SPatrick Sanan Fortran Notes: 99*95452b02SPatrick Sanan Only a single monitor function can be set for each SNESLineSearch object 100dcf2fd19SBarry Smith 101dcf2fd19SBarry Smith Level: intermediate 102dcf2fd19SBarry Smith 103dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 104dcf2fd19SBarry Smith 105dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel() 106dcf2fd19SBarry Smith @*/ 107dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 108dcf2fd19SBarry Smith { 10978064530SBarry Smith PetscErrorCode ierr; 11078064530SBarry Smith PetscInt i; 11178064530SBarry Smith PetscBool identical; 11278064530SBarry Smith 113dcf2fd19SBarry Smith PetscFunctionBegin; 114dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 11578064530SBarry Smith for (i=0; i<ls->numbermonitors;i++) { 11678064530SBarry Smith ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr); 11778064530SBarry Smith if (identical) PetscFunctionReturn(0); 11878064530SBarry Smith } 119dcf2fd19SBarry Smith if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 120dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 121dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 122dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void*)mctx; 123dcf2fd19SBarry Smith PetscFunctionReturn(0); 124dcf2fd19SBarry Smith } 125dcf2fd19SBarry Smith 126dcf2fd19SBarry Smith /*@C 127dcf2fd19SBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries 128dcf2fd19SBarry Smith 129dcf2fd19SBarry Smith Collective on SNESLineSearch 130dcf2fd19SBarry Smith 131dcf2fd19SBarry Smith Input Parameters: 132dcf2fd19SBarry Smith + ls - the SNES linesearch object 133d12e167eSBarry Smith - vf - the context for the monitor, in this case it is an ASCII PetscViewer and format 134dcf2fd19SBarry Smith 135dcf2fd19SBarry Smith Level: intermediate 136dcf2fd19SBarry Smith 137dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm 138dcf2fd19SBarry Smith 139dcf2fd19SBarry Smith .seealso: SNESMonitorSet(), SNESMonitorSolution() 140dcf2fd19SBarry Smith @*/ 141d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf) 142dcf2fd19SBarry Smith { 143dcf2fd19SBarry Smith PetscErrorCode ierr; 144d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 145dcf2fd19SBarry Smith Vec Y,W,G; 146dcf2fd19SBarry Smith 147dcf2fd19SBarry Smith PetscFunctionBegin; 148dcf2fd19SBarry Smith ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr); 149d12e167eSBarry Smith ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 150dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr); 151dcf2fd19SBarry Smith ierr = VecView(Y,viewer);CHKERRQ(ierr); 152dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr); 153dcf2fd19SBarry Smith ierr = VecView(W,viewer);CHKERRQ(ierr); 154dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr); 155dcf2fd19SBarry Smith ierr = VecView(G,viewer);CHKERRQ(ierr); 156d12e167eSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 157dcf2fd19SBarry Smith PetscFunctionReturn(0); 158dcf2fd19SBarry Smith } 159dcf2fd19SBarry Smith 160f40b411bSPeter Brune /*@ 161cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 162f40b411bSPeter Brune 163cd7522eaSPeter Brune Logically Collective on Comm 164f40b411bSPeter Brune 165f40b411bSPeter Brune Input Parameters: 166cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 167f40b411bSPeter Brune 168f40b411bSPeter Brune Output Parameters: 1698e557f58SPeter Brune . outlinesearch - the new linesearch context 170f40b411bSPeter Brune 171162e0bf5SPeter Brune Level: developer 172162e0bf5SPeter Brune 173162e0bf5SPeter Brune Notes: 174162e0bf5SPeter Brune The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance 175162e0bf5SPeter Brune already associated with the SNES. This function is for developer use. 176f40b411bSPeter Brune 177cd7522eaSPeter Brune .keywords: LineSearch, create, context 178f40b411bSPeter Brune 179162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch() 180f40b411bSPeter Brune @*/ 181f40b411bSPeter Brune 182bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 183bf388a1fSBarry Smith { 184bf7f4e0aSPeter Brune PetscErrorCode ierr; 185f1c6b773SPeter Brune SNESLineSearch linesearch; 186bf388a1fSBarry Smith 187bf7f4e0aSPeter Brune PetscFunctionBegin; 188ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 189f34a81feSMatthew G. Knepley ierr = SNESInitializePackage();CHKERRQ(ierr); 1900298fd71SBarry Smith *outlinesearch = NULL; 191f5af7f23SKarl Rupp 19273107ff1SLisandro Dalcin ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 193bf7f4e0aSPeter Brune 1940298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1950298fd71SBarry Smith linesearch->vec_func_new = NULL; 1960298fd71SBarry Smith linesearch->vec_sol = NULL; 1970298fd71SBarry Smith linesearch->vec_func = NULL; 1980298fd71SBarry Smith linesearch->vec_update = NULL; 1999bd66eb0SPeter Brune 200bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 201bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 202bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 203bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 204422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 205bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 206bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 207bf7f4e0aSPeter Brune linesearch->damping = 1.0; 208bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 209bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 210516fe3c3SPeter Brune linesearch->rtol = 1e-8; 211516fe3c3SPeter Brune linesearch->atol = 1e-15; 212516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2130298fd71SBarry Smith linesearch->precheckctx = NULL; 2140298fd71SBarry Smith linesearch->postcheckctx = NULL; 21559405d5eSPeter Brune linesearch->max_its = 1; 216bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 217bf7f4e0aSPeter Brune *outlinesearch = linesearch; 218bf7f4e0aSPeter Brune PetscFunctionReturn(0); 219bf7f4e0aSPeter Brune } 220bf7f4e0aSPeter Brune 221f40b411bSPeter Brune /*@ 22278bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 22378bcb3b5SPeter Brune any required vectors. 224f40b411bSPeter Brune 225cd7522eaSPeter Brune Collective on SNESLineSearch 226f40b411bSPeter Brune 227f40b411bSPeter Brune Input Parameters: 228f40b411bSPeter Brune . linesearch - The LineSearch instance. 229f40b411bSPeter Brune 230cd7522eaSPeter Brune Notes: 231f190f2fcSBarry Smith For most cases, this needn't be called by users or outside of SNESLineSearchApply(). 232cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 23378bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 234cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 235cd7522eaSPeter Brune allocated upfront. 236cd7522eaSPeter Brune 23778bcb3b5SPeter Brune Level: advanced 238f40b411bSPeter Brune 239f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 240f40b411bSPeter Brune 241f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 242f40b411bSPeter Brune @*/ 243f40b411bSPeter Brune 244bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 245bf388a1fSBarry Smith { 246bf7f4e0aSPeter Brune PetscErrorCode ierr; 247bf388a1fSBarry Smith 248bf7f4e0aSPeter Brune PetscFunctionBegin; 249bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 2501a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 251bf7f4e0aSPeter Brune } 252bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 2539bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 254bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 2559bd66eb0SPeter Brune } 2569bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 2579bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 2589bd66eb0SPeter Brune } 259bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 260bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 261bf7f4e0aSPeter Brune } 262ed07d7d7SPeter Brune if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);} 263bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 264bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 265bf7f4e0aSPeter Brune } 266bf7f4e0aSPeter Brune PetscFunctionReturn(0); 267bf7f4e0aSPeter Brune } 268bf7f4e0aSPeter Brune 269f40b411bSPeter Brune 270f40b411bSPeter Brune /*@ 271f190f2fcSBarry Smith SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search. 272f40b411bSPeter Brune 273f1c6b773SPeter Brune Collective on SNESLineSearch 274f40b411bSPeter Brune 275f40b411bSPeter Brune Input Parameters: 276f40b411bSPeter Brune . linesearch - The LineSearch instance. 277f40b411bSPeter Brune 278*95452b02SPatrick Sanan Notes: 279*95452b02SPatrick Sanan Usually only called by SNESReset() 280f190f2fcSBarry Smith 281f190f2fcSBarry Smith Level: developer 282f40b411bSPeter Brune 283cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 284f40b411bSPeter Brune 285f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 286f40b411bSPeter Brune @*/ 287f40b411bSPeter Brune 288bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 289bf388a1fSBarry Smith { 290bf7f4e0aSPeter Brune PetscErrorCode ierr; 291bf388a1fSBarry Smith 292bf7f4e0aSPeter Brune PetscFunctionBegin; 293f5af7f23SKarl Rupp if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch); 294f5af7f23SKarl Rupp 295bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 296bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 297bf7f4e0aSPeter Brune 298bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 299f5af7f23SKarl Rupp 300bf7f4e0aSPeter Brune linesearch->nwork = 0; 301bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 302bf7f4e0aSPeter Brune PetscFunctionReturn(0); 303bf7f4e0aSPeter Brune } 304bf7f4e0aSPeter Brune 305ed07d7d7SPeter Brune /*@C 306f190f2fcSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search 307ed07d7d7SPeter Brune 308ed07d7d7SPeter Brune Input Parameters: 309ed07d7d7SPeter Brune . linesearch - the SNESLineSearch context 310f190f2fcSBarry Smith + func - function evaluation routine 311ed07d7d7SPeter Brune 312ed07d7d7SPeter Brune Level: developer 313ed07d7d7SPeter Brune 314*95452b02SPatrick Sanan Notes: 315*95452b02SPatrick Sanan This is used internally by PETSc and not called by users 316f190f2fcSBarry Smith 317ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check 318ed07d7d7SPeter Brune 319f190f2fcSBarry Smith .seealso: SNESSetFunction() 320ed07d7d7SPeter Brune @*/ 321ed07d7d7SPeter Brune PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec)) 322ed07d7d7SPeter Brune { 323ed07d7d7SPeter Brune PetscFunctionBegin; 324ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 325ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 326ed07d7d7SPeter Brune PetscFunctionReturn(0); 327ed07d7d7SPeter Brune } 328ed07d7d7SPeter Brune 329ed07d7d7SPeter Brune 3306b2b7091SBarry Smith /*MC 331f190f2fcSBarry Smith SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called 3326b2b7091SBarry Smith 3336b2b7091SBarry Smith Synopsis: 334aaa7dc30SBarry Smith #include <petscsnes.h> 3356b2b7091SBarry Smith SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 3366b2b7091SBarry Smith 3376b2b7091SBarry Smith Input Parameters: 3386b2b7091SBarry Smith + x - solution vector 3396b2b7091SBarry Smith . y - search direction vector 3406b2b7091SBarry Smith - changed - flag to indicate the precheck changed x or y. 3416b2b7091SBarry Smith 342f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck() 343f190f2fcSBarry Smith and SNESLineSearchGetPreCheck() 344f190f2fcSBarry Smith 345878cb397SSatish Balay Level: advanced 346878cb397SSatish Balay 347f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck() 3486b2b7091SBarry Smith M*/ 3496b2b7091SBarry Smith 35086d74e61SPeter Brune /*@C 351f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 352df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 353f190f2fcSBarry Smith determined the search direction. 35486d74e61SPeter Brune 355f1c6b773SPeter Brune Logically Collective on SNESLineSearch 35686d74e61SPeter Brune 35786d74e61SPeter Brune Input Parameters: 358f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 359f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence 360f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 36186d74e61SPeter Brune 36286d74e61SPeter Brune Level: intermediate 36386d74e61SPeter Brune 36486d74e61SPeter Brune .keywords: set, linesearch, pre-check 36586d74e61SPeter Brune 366f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 36786d74e61SPeter Brune @*/ 368f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx) 36986d74e61SPeter Brune { 3709bd66eb0SPeter Brune PetscFunctionBegin; 371f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 372f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 37386d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 37486d74e61SPeter Brune PetscFunctionReturn(0); 37586d74e61SPeter Brune } 37686d74e61SPeter Brune 37786d74e61SPeter Brune /*@C 37853deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 37986d74e61SPeter Brune 38086d74e61SPeter Brune Input Parameters: 381f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 38286d74e61SPeter Brune 38386d74e61SPeter Brune Output Parameters: 384f190f2fcSBarry Smith + func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence 385f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 38686d74e61SPeter Brune 38786d74e61SPeter Brune Level: intermediate 38886d74e61SPeter Brune 38986d74e61SPeter Brune .keywords: get, linesearch, pre-check 39086d74e61SPeter Brune 391f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 39286d74e61SPeter Brune @*/ 393f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx) 39486d74e61SPeter Brune { 3959bd66eb0SPeter Brune PetscFunctionBegin; 396f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 397f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 39886d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 39986d74e61SPeter Brune PetscFunctionReturn(0); 40086d74e61SPeter Brune } 40186d74e61SPeter Brune 4026b2b7091SBarry Smith /*MC 403f190f2fcSBarry Smith SNESLineSearchPostCheckFunction - form of function that is called after line search is complete 4046b2b7091SBarry Smith 4056b2b7091SBarry Smith Synopsis: 406aaa7dc30SBarry Smith #include <petscsnes.h> 4076b2b7091SBarry Smith SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y, Vec w, *changed_y, PetscBool *changed_w); 4086b2b7091SBarry Smith 4096b2b7091SBarry Smith Input Parameters: 4106b2b7091SBarry Smith + x - old solution vector 4116b2b7091SBarry Smith . y - search direction vector 4126b2b7091SBarry Smith . w - new solution vector 4136b2b7091SBarry Smith . changed_y - indicates that the line search changed y 4146b2b7091SBarry Smith - changed_w - indicates that the line search changed w 4156b2b7091SBarry Smith 416f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck() 417f190f2fcSBarry Smith and SNESLineSearchGetPostCheck() 418f190f2fcSBarry Smith 419878cb397SSatish Balay Level: advanced 4206b2b7091SBarry Smith 421f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck() 4226b2b7091SBarry Smith M*/ 4236b2b7091SBarry Smith 42486d74e61SPeter Brune /*@C 425f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 426f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 42786d74e61SPeter Brune 428f1c6b773SPeter Brune Logically Collective on SNESLineSearch 42986d74e61SPeter Brune 43086d74e61SPeter Brune Input Parameters: 431f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 432f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence 433f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 43486d74e61SPeter Brune 43586d74e61SPeter Brune Level: intermediate 43686d74e61SPeter Brune 43786d74e61SPeter Brune .keywords: set, linesearch, post-check 43886d74e61SPeter Brune 439f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 44086d74e61SPeter Brune @*/ 441f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 44286d74e61SPeter Brune { 44386d74e61SPeter Brune PetscFunctionBegin; 444f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 445f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 44686d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 44786d74e61SPeter Brune PetscFunctionReturn(0); 44886d74e61SPeter Brune } 44986d74e61SPeter Brune 45086d74e61SPeter Brune /*@C 451f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 45286d74e61SPeter Brune 45386d74e61SPeter Brune Input Parameters: 454f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 45586d74e61SPeter Brune 45686d74e61SPeter Brune Output Parameters: 457f190f2fcSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction 458f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 45986d74e61SPeter Brune 46086d74e61SPeter Brune Level: intermediate 46186d74e61SPeter Brune 46286d74e61SPeter Brune .keywords: get, linesearch, post-check 46386d74e61SPeter Brune 464f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck() 46586d74e61SPeter Brune @*/ 466f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 46786d74e61SPeter Brune { 4689bd66eb0SPeter Brune PetscFunctionBegin; 469f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 470f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 47186d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 47286d74e61SPeter Brune PetscFunctionReturn(0); 47386d74e61SPeter Brune } 47486d74e61SPeter Brune 475f40b411bSPeter Brune /*@ 476f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 477f40b411bSPeter Brune 478cd7522eaSPeter Brune Logically Collective on SNESLineSearch 479f40b411bSPeter Brune 480f40b411bSPeter Brune Input Parameters: 4817b1df9c1SPeter Brune + linesearch - The linesearch instance. 4827b1df9c1SPeter Brune . X - The current solution 4837b1df9c1SPeter Brune - Y - The step direction 484f40b411bSPeter Brune 485f40b411bSPeter Brune Output Parameters: 4868e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 487f40b411bSPeter Brune 488f190f2fcSBarry Smith Level: developer 489f40b411bSPeter Brune 490f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 491f40b411bSPeter Brune 492f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 493f40b411bSPeter Brune @*/ 4947b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 495bf7f4e0aSPeter Brune { 496bf7f4e0aSPeter Brune PetscErrorCode ierr; 4975fd66863SKarl Rupp 498bf7f4e0aSPeter Brune PetscFunctionBegin; 499bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 5006b2b7091SBarry Smith if (linesearch->ops->precheck) { 5016b2b7091SBarry Smith ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr); 50238bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 503bf7f4e0aSPeter Brune } 504bf7f4e0aSPeter Brune PetscFunctionReturn(0); 505bf7f4e0aSPeter Brune } 506bf7f4e0aSPeter Brune 507f40b411bSPeter Brune /*@ 508f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 509f40b411bSPeter Brune 510cd7522eaSPeter Brune Logically Collective on SNESLineSearch 511f40b411bSPeter Brune 512f40b411bSPeter Brune Input Parameters: 5137b1df9c1SPeter Brune + linesearch - The linesearch context 5147b1df9c1SPeter Brune . X - The last solution 5157b1df9c1SPeter Brune . Y - The step direction 5167b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 517f40b411bSPeter Brune 518f40b411bSPeter Brune Output Parameters: 51978bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 52078bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 521f40b411bSPeter Brune 522f190f2fcSBarry Smith Level: developer 523f40b411bSPeter Brune 524f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 525f40b411bSPeter Brune 526f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 527f40b411bSPeter Brune @*/ 5287b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 529bf7f4e0aSPeter Brune { 530bf7f4e0aSPeter Brune PetscErrorCode ierr; 531bf388a1fSBarry Smith 532bf7f4e0aSPeter Brune PetscFunctionBegin; 533bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 534bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 5356b2b7091SBarry Smith if (linesearch->ops->postcheck) { 5366b2b7091SBarry Smith ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr); 53738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 53838bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 53986d74e61SPeter Brune } 54086d74e61SPeter Brune PetscFunctionReturn(0); 54186d74e61SPeter Brune } 54286d74e61SPeter Brune 54386d74e61SPeter Brune /*@C 54486d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 54586d74e61SPeter Brune 546cd7522eaSPeter Brune Logically Collective on SNESLineSearch 54786d74e61SPeter Brune 54886d74e61SPeter Brune Input Arguments: 54986d74e61SPeter Brune + linesearch - linesearch context 55086d74e61SPeter Brune . X - base state for this step 55186d74e61SPeter Brune . Y - initial correction 552907376e6SBarry Smith - ctx - context for this function 55386d74e61SPeter Brune 55486d74e61SPeter Brune Output Arguments: 55586d74e61SPeter Brune + Y - correction, possibly modified 55686d74e61SPeter Brune - changed - flag indicating that Y was modified 55786d74e61SPeter Brune 55886d74e61SPeter Brune Options Database Key: 559cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 560cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 56186d74e61SPeter Brune 56286d74e61SPeter Brune Level: advanced 56386d74e61SPeter Brune 56486d74e61SPeter Brune Notes: 56586d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 56686d74e61SPeter Brune 56786d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 56886d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 56986d74e61SPeter Brune is generally not useful when using a Newton linearization. 57086d74e61SPeter Brune 57186d74e61SPeter Brune Reference: 57286d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 57386d74e61SPeter Brune 57486d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 57586d74e61SPeter Brune @*/ 576f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 57786d74e61SPeter Brune { 57886d74e61SPeter Brune PetscErrorCode ierr; 57986d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 58086d74e61SPeter Brune Vec Ylast; 58186d74e61SPeter Brune PetscScalar dot; 58286d74e61SPeter Brune PetscInt iter; 58386d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 58486d74e61SPeter Brune SNES snes; 58586d74e61SPeter Brune 58686d74e61SPeter Brune PetscFunctionBegin; 587f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 58886d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 58986d74e61SPeter Brune if (!Ylast) { 59086d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 59186d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 59286d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 59386d74e61SPeter Brune } 59486d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 59586d74e61SPeter Brune if (iter < 2) { 59686d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 59786d74e61SPeter Brune *changed = PETSC_FALSE; 59886d74e61SPeter Brune PetscFunctionReturn(0); 59986d74e61SPeter Brune } 60086d74e61SPeter Brune 60186d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 60286d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 60386d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 60486d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 605255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 60686d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 60786d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 60886d74e61SPeter Brune /* Modify the step Y */ 60986d74e61SPeter Brune PetscReal alpha,ydiffnorm; 61086d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 61186d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 61286d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 61386d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 61486d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 615c69d1a72SBarry 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); 61686d74e61SPeter Brune } else { 617c69d1a72SBarry 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); 61886d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 61986d74e61SPeter Brune *changed = PETSC_FALSE; 620bf7f4e0aSPeter Brune } 621bf7f4e0aSPeter Brune PetscFunctionReturn(0); 622bf7f4e0aSPeter Brune } 623bf7f4e0aSPeter Brune 624f40b411bSPeter Brune /*@ 625cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 626f40b411bSPeter Brune 627f1c6b773SPeter Brune Collective on SNESLineSearch 628f40b411bSPeter Brune 629f40b411bSPeter Brune Input Parameters: 6308e557f58SPeter Brune + linesearch - The linesearch context 6318e557f58SPeter Brune . X - The current solution 6328e557f58SPeter Brune . F - The current function 6338e557f58SPeter Brune . fnorm - The current norm 6348e557f58SPeter Brune - Y - The search direction 635f40b411bSPeter Brune 636f40b411bSPeter Brune Output Parameters: 6378e557f58SPeter Brune + X - The new solution 6388e557f58SPeter Brune . F - The new function 6398e557f58SPeter Brune - fnorm - The new function norm 640f40b411bSPeter Brune 641cd7522eaSPeter Brune Options Database Keys: 642d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell 643dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 6441fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 6451fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 6463c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 6473c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 648cd7522eaSPeter Brune 649cd7522eaSPeter Brune Notes: 650cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 651cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 652cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 653cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 654cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 655cd7522eaSPeter Brune therefore may be fairly expensive. 656cd7522eaSPeter Brune 657f40b411bSPeter Brune Level: Intermediate 658f40b411bSPeter Brune 659f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 660f40b411bSPeter Brune 6611fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(), 6621fe24845SBarry Smith SNESLineSearchType, SNESLineSearchSetType() 663f40b411bSPeter Brune @*/ 664bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 665bf388a1fSBarry Smith { 666bf7f4e0aSPeter Brune PetscErrorCode ierr; 667bf7f4e0aSPeter Brune 668bf388a1fSBarry Smith PetscFunctionBegin; 669f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 670bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 671bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 672bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 673bf7f4e0aSPeter Brune 674422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 675bf7f4e0aSPeter Brune 676bf7f4e0aSPeter Brune linesearch->vec_sol = X; 677bf7f4e0aSPeter Brune linesearch->vec_update = Y; 678bf7f4e0aSPeter Brune linesearch->vec_func = F; 679bf7f4e0aSPeter Brune 680f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 681bf7f4e0aSPeter Brune 682f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 683bf7f4e0aSPeter Brune 684f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 685f5af7f23SKarl Rupp else { 686bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 687bf7f4e0aSPeter Brune } 688bf7f4e0aSPeter Brune 68957a83faaSBarry Smith ierr = PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 690bf7f4e0aSPeter Brune 691bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 692bf7f4e0aSPeter Brune 69357a83faaSBarry Smith ierr = PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 694bf7f4e0aSPeter Brune 695f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 696bf7f4e0aSPeter Brune PetscFunctionReturn(0); 697bf7f4e0aSPeter Brune } 698bf7f4e0aSPeter Brune 699f40b411bSPeter Brune /*@ 700f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 701f40b411bSPeter Brune 702f1c6b773SPeter Brune Collective on SNESLineSearch 703f40b411bSPeter Brune 704f40b411bSPeter Brune Input Parameters: 7058e557f58SPeter Brune . linesearch - The linesearch context 706f40b411bSPeter Brune 707f40b411bSPeter Brune Level: Intermediate 708f40b411bSPeter Brune 70978bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 710f40b411bSPeter Brune 711cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 712f40b411bSPeter Brune @*/ 713bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 714bf388a1fSBarry Smith { 715bf7f4e0aSPeter Brune PetscErrorCode ierr; 716bf388a1fSBarry Smith 717bf7f4e0aSPeter Brune PetscFunctionBegin; 718bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 719f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 720bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 721e04113cfSBarry Smith ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr); 72222d28d08SBarry Smith ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr); 723f5af7f23SKarl Rupp if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch); 724bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 725dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr); 726e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 727bf7f4e0aSPeter Brune PetscFunctionReturn(0); 728bf7f4e0aSPeter Brune } 729bf7f4e0aSPeter Brune 730f40b411bSPeter Brune /*@ 731dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 732bf7f4e0aSPeter Brune 733bf7f4e0aSPeter Brune Input Parameters: 734dcf2fd19SBarry Smith + linesearch - the linesearch object 735dcf2fd19SBarry Smith - viewer - an ASCII PetscViewer or NULL to turn off monitor 736bf7f4e0aSPeter Brune 737dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 738bf7f4e0aSPeter Brune 739bf7f4e0aSPeter Brune Options Database: 740dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 741bf7f4e0aSPeter Brune 742bf7f4e0aSPeter Brune Level: intermediate 743bf7f4e0aSPeter Brune 744dcf2fd19SBarry Smith Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with 745d12e167eSBarry Smith SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the 746d12e167eSBarry Smith line search that are not visible to the other monitors. 747dcf2fd19SBarry Smith 748d12e167eSBarry Smith .seealso: SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor() 749bf7f4e0aSPeter Brune @*/ 750dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 751bf7f4e0aSPeter Brune { 752bf7f4e0aSPeter Brune PetscErrorCode ierr; 753bf388a1fSBarry Smith 754bf7f4e0aSPeter Brune PetscFunctionBegin; 755dcf2fd19SBarry Smith if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);} 756bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 757dcf2fd19SBarry Smith linesearch->monitor = viewer; 758bf7f4e0aSPeter Brune PetscFunctionReturn(0); 759bf7f4e0aSPeter Brune } 760bf7f4e0aSPeter Brune 761f40b411bSPeter Brune /*@ 762dcf2fd19SBarry Smith SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor. 7636a388c36SPeter Brune 764f190f2fcSBarry Smith Input Parameter: 7658e557f58SPeter Brune . linesearch - linesearch context 766f40b411bSPeter Brune 767f190f2fcSBarry Smith Output Parameter: 7688e557f58SPeter Brune . monitor - monitor context 769f40b411bSPeter Brune 770f40b411bSPeter Brune Logically Collective on SNES 771f40b411bSPeter Brune 772f40b411bSPeter Brune Options Database Keys: 7738e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 774f40b411bSPeter Brune 775f40b411bSPeter Brune Level: intermediate 776f40b411bSPeter Brune 777dcf2fd19SBarry Smith .seealso: SNESLineSearchSetDefaultMonitor(), PetscViewer 778f40b411bSPeter Brune @*/ 779dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 7806a388c36SPeter Brune { 7816a388c36SPeter Brune PetscFunctionBegin; 782f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 7836a388c36SPeter Brune if (monitor) { 7846a388c36SPeter Brune PetscValidPointer(monitor, 2); 7856a388c36SPeter Brune *monitor = linesearch->monitor; 7866a388c36SPeter Brune } 7876a388c36SPeter Brune PetscFunctionReturn(0); 7886a388c36SPeter Brune } 7896a388c36SPeter Brune 790dcf2fd19SBarry Smith /*@C 791dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 792dcf2fd19SBarry Smith 793dcf2fd19SBarry Smith Collective on SNESLineSearch 794dcf2fd19SBarry Smith 795dcf2fd19SBarry Smith Input Parameters: 796dcf2fd19SBarry Smith + ls - LineSearch object you wish to monitor 797dcf2fd19SBarry Smith . name - the monitor type one is seeking 798dcf2fd19SBarry Smith . help - message indicating what monitoring is done 799dcf2fd19SBarry Smith . manual - manual page for the monitor 800dcf2fd19SBarry Smith . monitor - the monitor function 801dcf2fd19SBarry 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 802dcf2fd19SBarry Smith 803dcf2fd19SBarry Smith Level: developer 804dcf2fd19SBarry Smith 805dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 806dcf2fd19SBarry Smith PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 807dcf2fd19SBarry Smith PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 808dcf2fd19SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 809dcf2fd19SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 810dcf2fd19SBarry Smith PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 811dcf2fd19SBarry Smith PetscOptionsFList(), PetscOptionsEList() 812dcf2fd19SBarry Smith @*/ 813d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*)) 814dcf2fd19SBarry Smith { 815dcf2fd19SBarry Smith PetscErrorCode ierr; 816dcf2fd19SBarry Smith PetscViewer viewer; 817dcf2fd19SBarry Smith PetscViewerFormat format; 818dcf2fd19SBarry Smith PetscBool flg; 819dcf2fd19SBarry Smith 820dcf2fd19SBarry Smith PetscFunctionBegin; 821dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 822dcf2fd19SBarry Smith if (flg) { 823d12e167eSBarry Smith PetscViewerAndFormat *vf; 824d12e167eSBarry Smith ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 825d12e167eSBarry Smith ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 826dcf2fd19SBarry Smith if (monitorsetup) { 827d12e167eSBarry Smith ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr); 828dcf2fd19SBarry Smith } 829d12e167eSBarry Smith ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 830dcf2fd19SBarry Smith } 831dcf2fd19SBarry Smith PetscFunctionReturn(0); 832dcf2fd19SBarry Smith } 833dcf2fd19SBarry Smith 834f40b411bSPeter Brune /*@ 835f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 836f40b411bSPeter Brune 837f40b411bSPeter Brune Input Parameters: 8388e557f58SPeter Brune . linesearch - linesearch context 839f40b411bSPeter Brune 840f40b411bSPeter Brune Options Database Keys: 841d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell 8423c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 8431fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms()) 84471eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 8451a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 8461a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 8471a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 8481a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 8491a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 850dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 851dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 8528e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 853cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 8541a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 8551a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method 856f40b411bSPeter Brune 857f1c6b773SPeter Brune Logically Collective on SNESLineSearch 858f40b411bSPeter Brune 859f40b411bSPeter Brune Level: intermediate 860f40b411bSPeter Brune 8611fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(), 8621fe24845SBarry Smith SNESLineSearchType, SNESLineSearchSetComputeNorms() 863f40b411bSPeter Brune @*/ 864bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 865bf388a1fSBarry Smith { 866bf7f4e0aSPeter Brune PetscErrorCode ierr; 8671a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 868bf7f4e0aSPeter Brune char type[256]; 869bf7f4e0aSPeter Brune PetscBool flg, set; 870dcf2fd19SBarry Smith PetscViewer viewer; 871bf388a1fSBarry Smith 872bf7f4e0aSPeter Brune PetscFunctionBegin; 8730f51fdf8SToby Isaac ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr); 874bf7f4e0aSPeter Brune 875bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 876f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 877a264d7a6SBarry Smith ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 878bf7f4e0aSPeter Brune if (flg) { 879f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 880bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 881f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 882bf7f4e0aSPeter Brune } 883bf7f4e0aSPeter Brune 884dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr); 885dcf2fd19SBarry Smith if (set) { 886dcf2fd19SBarry Smith ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr); 887dcf2fd19SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 888dcf2fd19SBarry Smith } 889dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr); 890bf7f4e0aSPeter Brune 8911a9b3a06SPeter Brune /* tolerances */ 89294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr); 89394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr); 89494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr); 89594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr); 89694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr); 89794ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr); 8987a35526eSPeter Brune 8991a9b3a06SPeter Brune /* damping parameters */ 90094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr); 9011a9b3a06SPeter Brune 90294ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr); 9031a9b3a06SPeter Brune 9041a9b3a06SPeter Brune /* precheck */ 9057a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 9067a35526eSPeter Brune if (set) { 9077a35526eSPeter Brune if (flg) { 9087a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 909f5af7f23SKarl Rupp 9107a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 9110298fd71SBarry Smith "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr); 9127a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 9137a35526eSPeter Brune } else { 9140298fd71SBarry Smith ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr); 9157a35526eSPeter Brune } 9167a35526eSPeter Brune } 91794ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr); 91894ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr); 9197a35526eSPeter Brune 9205a9b6599SPeter Brune if (linesearch->ops->setfromoptions) { 921e55864a3SBarry Smith (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr); 9225a9b6599SPeter Brune } 9235a9b6599SPeter Brune 9240633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr); 925bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 926bf7f4e0aSPeter Brune PetscFunctionReturn(0); 927bf7f4e0aSPeter Brune } 928bf7f4e0aSPeter Brune 929f40b411bSPeter Brune /*@ 930f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 931f40b411bSPeter Brune 932f40b411bSPeter Brune Input Parameters: 9338e557f58SPeter Brune . linesearch - linesearch context 934f40b411bSPeter Brune 935f1c6b773SPeter Brune Logically Collective on SNESLineSearch 936f40b411bSPeter Brune 937f40b411bSPeter Brune Level: intermediate 938f40b411bSPeter Brune 939dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate() 940f40b411bSPeter Brune @*/ 941bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 942bf388a1fSBarry Smith { 9437f1410a3SPeter Brune PetscErrorCode ierr; 9447f1410a3SPeter Brune PetscBool iascii; 945bf388a1fSBarry Smith 946bf7f4e0aSPeter Brune PetscFunctionBegin; 9477f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9487f1410a3SPeter Brune if (!viewer) { 949ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr); 9507f1410a3SPeter Brune } 9517f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 9527f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 953f40b411bSPeter Brune 954251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9557f1410a3SPeter Brune if (iascii) { 95663b15415SAlp Dener PetscInt tabs; 95763b15415SAlp Dener ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr); 95863b15415SAlp Dener ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 959dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr); 9607f1410a3SPeter Brune if (linesearch->ops->view) { 9617f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9627f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 9637f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 9647f1410a3SPeter Brune } 965c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr); 966c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr); 9677f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 9686b2b7091SBarry Smith if (linesearch->ops->precheck) { 9696b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 9707f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 9717f1410a3SPeter Brune } else { 9727f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 9737f1410a3SPeter Brune } 9747f1410a3SPeter Brune } 9756b2b7091SBarry Smith if (linesearch->ops->postcheck) { 9767f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 9777f1410a3SPeter Brune } 97863b15415SAlp Dener ierr = PetscViewerASCIISetTab(viewer, tabs);CHKERRQ(ierr); 9797f1410a3SPeter Brune } 980bf7f4e0aSPeter Brune PetscFunctionReturn(0); 981bf7f4e0aSPeter Brune } 982bf7f4e0aSPeter Brune 983ea5d4fccSPeter Brune /*@C 984f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 985f40b411bSPeter Brune 986f190f2fcSBarry Smith Logically Collective on SNESLineSearch 987f190f2fcSBarry Smith 988f40b411bSPeter Brune Input Parameters: 9898e557f58SPeter Brune + linesearch - linesearch context 990f40b411bSPeter Brune - type - The type of line search to be used 991f40b411bSPeter Brune 992cd7522eaSPeter Brune Available Types: 9931fe24845SBarry Smith + SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step 9941fe24845SBarry Smith . SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function 9951fe24845SBarry Smith . SNESLINESEARCHL2 - Secant line search over the L2 norm of the function 9961fe24845SBarry Smith . SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 9971fe24845SBarry Smith . SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch 9981fe24845SBarry Smith - SNESLINESEARCHSHELL - User provided SNESLineSearch implementation 9991fe24845SBarry Smith 10001fe24845SBarry Smith Options Database: 10011fe24845SBarry Smith . -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell 1002cd7522eaSPeter Brune 1003f40b411bSPeter Brune Level: intermediate 1004f40b411bSPeter Brune 10051fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions() 1006f40b411bSPeter Brune @*/ 100719fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 1008bf7f4e0aSPeter Brune { 1009f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 1010bf7f4e0aSPeter Brune PetscBool match; 1011bf7f4e0aSPeter Brune 1012bf7f4e0aSPeter Brune PetscFunctionBegin; 1013f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1014bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 1015bf7f4e0aSPeter Brune 1016251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 1017bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 1018bf7f4e0aSPeter Brune 10191c9cd337SJed Brown ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr); 1020bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 1021bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 1022bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 1023bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 1024f5af7f23SKarl Rupp 10250298fd71SBarry Smith linesearch->ops->destroy = NULL; 1026bf7f4e0aSPeter Brune } 1027f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 1028bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 1029bf7f4e0aSPeter Brune linesearch->ops->view = 0; 1030bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 1031bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 1032bf7f4e0aSPeter Brune 1033bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 1034bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 1035bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1036bf7f4e0aSPeter Brune } 1037bf7f4e0aSPeter Brune 1038f40b411bSPeter Brune /*@ 103978bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 1040f40b411bSPeter Brune 1041f40b411bSPeter Brune Input Parameters: 10428e557f58SPeter Brune + linesearch - linesearch context 1043f40b411bSPeter Brune - snes - The snes instance 1044f40b411bSPeter Brune 104578bcb3b5SPeter Brune Level: developer 104678bcb3b5SPeter Brune 104778bcb3b5SPeter Brune Notes: 1048f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 10497601faf0SJed Brown SNESGetLineSearch(). This routine is therefore mainly called within SNES 105078bcb3b5SPeter Brune implementations. 1051f40b411bSPeter Brune 10528141a3b9SPeter Brune Level: developer 1053f40b411bSPeter Brune 1054cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1055f40b411bSPeter Brune @*/ 1056bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1057bf388a1fSBarry Smith { 1058bf7f4e0aSPeter Brune PetscFunctionBegin; 1059f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1060bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 1061bf7f4e0aSPeter Brune linesearch->snes = snes; 1062bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1063bf7f4e0aSPeter Brune } 1064bf7f4e0aSPeter Brune 1065f40b411bSPeter Brune /*@ 10668141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 10678141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 10688141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 10698141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 1070f40b411bSPeter Brune 1071f40b411bSPeter Brune Input Parameters: 10728e557f58SPeter Brune . linesearch - linesearch context 1073f40b411bSPeter Brune 1074f40b411bSPeter Brune Output Parameters: 1075f40b411bSPeter Brune . snes - The snes instance 1076f40b411bSPeter Brune 10778141a3b9SPeter Brune Level: developer 1078f40b411bSPeter Brune 1079cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1080f40b411bSPeter Brune @*/ 1081bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1082bf388a1fSBarry Smith { 1083bf7f4e0aSPeter Brune PetscFunctionBegin; 1084f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 10856a388c36SPeter Brune PetscValidPointer(snes, 2); 1086bf7f4e0aSPeter Brune *snes = linesearch->snes; 1087bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1088bf7f4e0aSPeter Brune } 1089bf7f4e0aSPeter Brune 1090f40b411bSPeter Brune /*@ 1091f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1092f40b411bSPeter Brune 1093f40b411bSPeter Brune Input Parameters: 10948e557f58SPeter Brune . linesearch - linesearch context 1095f40b411bSPeter Brune 1096f40b411bSPeter Brune Output Parameters: 1097cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 1098f40b411bSPeter Brune 109978bcb3b5SPeter Brune Level: advanced 110078bcb3b5SPeter Brune 11018e557f58SPeter Brune Notes: 11028e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 110378bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 110478bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 110578bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 110678bcb3b5SPeter Brune 1107cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 1108f40b411bSPeter Brune @*/ 1109f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 11106a388c36SPeter Brune { 11116a388c36SPeter Brune PetscFunctionBegin; 1112f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11136a388c36SPeter Brune PetscValidPointer(lambda, 2); 11146a388c36SPeter Brune *lambda = linesearch->lambda; 11156a388c36SPeter Brune PetscFunctionReturn(0); 11166a388c36SPeter Brune } 11176a388c36SPeter Brune 1118f40b411bSPeter Brune /*@ 1119f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 1120f40b411bSPeter Brune 1121f40b411bSPeter Brune Input Parameters: 11228e557f58SPeter Brune + linesearch - linesearch context 1123f40b411bSPeter Brune - lambda - The last steplength. 1124f40b411bSPeter Brune 1125cd7522eaSPeter Brune Notes: 1126f190f2fcSBarry Smith This routine is typically used within implementations of SNESLineSearchApply() 1127cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 1128cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1129cd7522eaSPeter Brune as an inner scaling parameter. 1130cd7522eaSPeter Brune 113178bcb3b5SPeter Brune Level: advanced 1132f40b411bSPeter Brune 1133f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 1134f40b411bSPeter Brune @*/ 1135f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 11366a388c36SPeter Brune { 11376a388c36SPeter Brune PetscFunctionBegin; 1138f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11396a388c36SPeter Brune linesearch->lambda = lambda; 11406a388c36SPeter Brune PetscFunctionReturn(0); 11416a388c36SPeter Brune } 11426a388c36SPeter Brune 1143f40b411bSPeter Brune /*@ 11443c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 114578bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 114678bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 114778bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1148f40b411bSPeter Brune 1149f40b411bSPeter Brune Input Parameters: 11508e557f58SPeter Brune . linesearch - linesearch context 1151f40b411bSPeter Brune 1152f40b411bSPeter Brune Output Parameters: 1153516fe3c3SPeter Brune + steptol - The minimum steplength 11546cc8e53bSPeter Brune . maxstep - The maximum steplength 1155516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1156516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1157516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1158516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1159f40b411bSPeter Brune 116078bcb3b5SPeter Brune Level: intermediate 116178bcb3b5SPeter Brune 116278bcb3b5SPeter Brune Notes: 116378bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 11643c7d6663SPeter Brune the type requires. 1165516fe3c3SPeter Brune 1166f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 1167f40b411bSPeter Brune @*/ 1168f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 11696a388c36SPeter Brune { 11706a388c36SPeter Brune PetscFunctionBegin; 1171f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1172516fe3c3SPeter Brune if (steptol) { 11736a388c36SPeter Brune PetscValidPointer(steptol, 2); 11746a388c36SPeter Brune *steptol = linesearch->steptol; 1175516fe3c3SPeter Brune } 1176516fe3c3SPeter Brune if (maxstep) { 1177516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 1178516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1179516fe3c3SPeter Brune } 1180516fe3c3SPeter Brune if (rtol) { 1181516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 1182516fe3c3SPeter Brune *rtol = linesearch->rtol; 1183516fe3c3SPeter Brune } 1184516fe3c3SPeter Brune if (atol) { 1185516fe3c3SPeter Brune PetscValidPointer(atol, 5); 1186516fe3c3SPeter Brune *atol = linesearch->atol; 1187516fe3c3SPeter Brune } 1188516fe3c3SPeter Brune if (ltol) { 1189516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 1190516fe3c3SPeter Brune *ltol = linesearch->ltol; 1191516fe3c3SPeter Brune } 1192516fe3c3SPeter Brune if (max_its) { 1193516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 1194516fe3c3SPeter Brune *max_its = linesearch->max_its; 1195516fe3c3SPeter Brune } 11966a388c36SPeter Brune PetscFunctionReturn(0); 11976a388c36SPeter Brune } 11986a388c36SPeter Brune 1199f40b411bSPeter Brune /*@ 12003c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 120178bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 120278bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 120378bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1204f40b411bSPeter Brune 1205f40b411bSPeter Brune Input Parameters: 12068e557f58SPeter Brune + linesearch - linesearch context 1207516fe3c3SPeter Brune . steptol - The minimum steplength 12086cc8e53bSPeter Brune . maxstep - The maximum steplength 1209516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1210516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1211516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1212516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1213f40b411bSPeter Brune 121478bcb3b5SPeter Brune Notes: 12153c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 121678bcb3b5SPeter Brune place of an argument. 1217f40b411bSPeter Brune 121878bcb3b5SPeter Brune Level: intermediate 1219516fe3c3SPeter Brune 1220f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 1221f40b411bSPeter Brune @*/ 1222f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 12236a388c36SPeter Brune { 12246a388c36SPeter Brune PetscFunctionBegin; 1225f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1226d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1227d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1228d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1229d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1230d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1231d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1232d3952184SSatish Balay 1233d3952184SSatish Balay if (steptol!= PETSC_DEFAULT) { 1234ce94432eSBarry Smith if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 12356a388c36SPeter Brune linesearch->steptol = steptol; 1236d3952184SSatish Balay } 1237d3952184SSatish Balay 1238d3952184SSatish Balay if (maxstep!= PETSC_DEFAULT) { 1239ce94432eSBarry Smith if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1240516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1241d3952184SSatish Balay } 1242d3952184SSatish Balay 1243d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1244ce94432eSBarry 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); 1245516fe3c3SPeter Brune linesearch->rtol = rtol; 1246d3952184SSatish Balay } 1247d3952184SSatish Balay 1248d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1249ce94432eSBarry Smith if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1250516fe3c3SPeter Brune linesearch->atol = atol; 1251d3952184SSatish Balay } 1252d3952184SSatish Balay 1253d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1254ce94432eSBarry Smith if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1255516fe3c3SPeter Brune linesearch->ltol = ltol; 1256d3952184SSatish Balay } 1257d3952184SSatish Balay 1258d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1259ce94432eSBarry Smith if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1260516fe3c3SPeter Brune linesearch->max_its = max_its; 1261d3952184SSatish Balay } 12626a388c36SPeter Brune PetscFunctionReturn(0); 12636a388c36SPeter Brune } 12646a388c36SPeter Brune 1265f40b411bSPeter Brune /*@ 1266f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1267f40b411bSPeter Brune 1268f40b411bSPeter Brune Input Parameters: 12698e557f58SPeter Brune . linesearch - linesearch context 1270f40b411bSPeter Brune 1271f40b411bSPeter Brune Output Parameters: 12728e557f58SPeter Brune . damping - The damping parameter 1273f40b411bSPeter Brune 127478bcb3b5SPeter Brune Level: advanced 1275f40b411bSPeter Brune 127678bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1277f40b411bSPeter Brune @*/ 1278f40b411bSPeter Brune 1279f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 12806a388c36SPeter Brune { 12816a388c36SPeter Brune PetscFunctionBegin; 1282f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12836a388c36SPeter Brune PetscValidPointer(damping, 2); 12846a388c36SPeter Brune *damping = linesearch->damping; 12856a388c36SPeter Brune PetscFunctionReturn(0); 12866a388c36SPeter Brune } 12876a388c36SPeter Brune 1288f40b411bSPeter Brune /*@ 1289f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1290f40b411bSPeter Brune 1291f40b411bSPeter Brune Input Parameters: 129203fd4120SBarry Smith + linesearch - linesearch context 129303fd4120SBarry Smith - damping - The damping parameter 1294f40b411bSPeter Brune 129503fd4120SBarry Smith Options Database: 129603fd4120SBarry Smith . -snes_linesearch_damping 1297f40b411bSPeter Brune Level: intermediate 1298f40b411bSPeter Brune 1299cd7522eaSPeter Brune Notes: 1300cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1301cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 130278bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1303cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1304cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1305cd7522eaSPeter Brune 1306f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1307f40b411bSPeter Brune @*/ 1308f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 13096a388c36SPeter Brune { 13106a388c36SPeter Brune PetscFunctionBegin; 1311f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13126a388c36SPeter Brune linesearch->damping = damping; 13136a388c36SPeter Brune PetscFunctionReturn(0); 13146a388c36SPeter Brune } 13156a388c36SPeter Brune 131659405d5eSPeter Brune /*@ 131759405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 131859405d5eSPeter Brune 131959405d5eSPeter Brune Input Parameters: 132078bcb3b5SPeter Brune . linesearch - linesearch context 132159405d5eSPeter Brune 132259405d5eSPeter Brune Output Parameters: 13238e557f58SPeter Brune . order - The order 132459405d5eSPeter Brune 132578bcb3b5SPeter Brune Possible Values for order: 13263c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 13273c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 13283c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 132978bcb3b5SPeter Brune 133059405d5eSPeter Brune Level: intermediate 133159405d5eSPeter Brune 133259405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 133359405d5eSPeter Brune @*/ 133459405d5eSPeter Brune 1335b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 133659405d5eSPeter Brune { 133759405d5eSPeter Brune PetscFunctionBegin; 133859405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 133959405d5eSPeter Brune PetscValidPointer(order, 2); 134059405d5eSPeter Brune *order = linesearch->order; 134159405d5eSPeter Brune PetscFunctionReturn(0); 134259405d5eSPeter Brune } 134359405d5eSPeter Brune 134459405d5eSPeter Brune /*@ 13451f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 134659405d5eSPeter Brune 134759405d5eSPeter Brune Input Parameters: 134878bcb3b5SPeter Brune . linesearch - linesearch context 134978bcb3b5SPeter Brune . order - The damping parameter 135059405d5eSPeter Brune 135159405d5eSPeter Brune Level: intermediate 135259405d5eSPeter Brune 135378bcb3b5SPeter Brune Possible Values for order: 13543c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 13553c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 13563c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 135778bcb3b5SPeter Brune 135859405d5eSPeter Brune Notes: 135959405d5eSPeter Brune Variable orders are supported by the following line searches: 136078bcb3b5SPeter Brune + bt - cubic and quadratic 136178bcb3b5SPeter Brune - cp - linear and quadratic 136259405d5eSPeter Brune 13631f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping() 136459405d5eSPeter Brune @*/ 1365b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 136659405d5eSPeter Brune { 136759405d5eSPeter Brune PetscFunctionBegin; 136859405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 136959405d5eSPeter Brune linesearch->order = order; 137059405d5eSPeter Brune PetscFunctionReturn(0); 137159405d5eSPeter Brune } 137259405d5eSPeter Brune 1373f40b411bSPeter Brune /*@ 1374f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1375f40b411bSPeter Brune 1376f40b411bSPeter Brune Input Parameters: 137778bcb3b5SPeter Brune . linesearch - linesearch context 1378f40b411bSPeter Brune 1379f40b411bSPeter Brune Output Parameters: 1380f40b411bSPeter Brune + xnorm - The norm of the current solution 1381f40b411bSPeter Brune . fnorm - The norm of the current function 1382f40b411bSPeter Brune - ynorm - The norm of the current update 1383f40b411bSPeter Brune 1384cd7522eaSPeter Brune Notes: 1385cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1386cd7522eaSPeter Brune 138778bcb3b5SPeter Brune Level: developer 1388f40b411bSPeter Brune 1389f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1390f40b411bSPeter Brune @*/ 1391f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1392bf7f4e0aSPeter Brune { 1393bf7f4e0aSPeter Brune PetscFunctionBegin; 1394f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1395f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1396f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1397f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 1398bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1399bf7f4e0aSPeter Brune } 1400bf7f4e0aSPeter Brune 1401f40b411bSPeter Brune /*@ 1402f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1403f40b411bSPeter Brune 1404f40b411bSPeter Brune Input Parameters: 140578bcb3b5SPeter Brune + linesearch - linesearch context 1406f40b411bSPeter Brune . xnorm - The norm of the current solution 1407f40b411bSPeter Brune . fnorm - The norm of the current function 1408f40b411bSPeter Brune - ynorm - The norm of the current update 1409f40b411bSPeter Brune 141078bcb3b5SPeter Brune Level: advanced 1411f40b411bSPeter Brune 1412f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1413f40b411bSPeter Brune @*/ 1414f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 14156a388c36SPeter Brune { 14166a388c36SPeter Brune PetscFunctionBegin; 1417f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14186a388c36SPeter Brune linesearch->xnorm = xnorm; 14196a388c36SPeter Brune linesearch->fnorm = fnorm; 14206a388c36SPeter Brune linesearch->ynorm = ynorm; 14216a388c36SPeter Brune PetscFunctionReturn(0); 14226a388c36SPeter Brune } 14236a388c36SPeter Brune 1424f40b411bSPeter Brune /*@ 1425f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1426f40b411bSPeter Brune 1427f40b411bSPeter Brune Input Parameters: 142878bcb3b5SPeter Brune . linesearch - linesearch context 1429f40b411bSPeter Brune 1430f40b411bSPeter Brune Options Database Keys: 14318e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1432f40b411bSPeter Brune 1433f40b411bSPeter Brune Level: intermediate 1434f40b411bSPeter Brune 143578bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1436f40b411bSPeter Brune @*/ 1437f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 14386a388c36SPeter Brune { 14396a388c36SPeter Brune PetscErrorCode ierr; 14409bd66eb0SPeter Brune SNES snes; 1441bf388a1fSBarry Smith 14426a388c36SPeter Brune PetscFunctionBegin; 14436a388c36SPeter Brune if (linesearch->norms) { 14449bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1445f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 14469bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 14479bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 14489bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 14499bd66eb0SPeter Brune } else { 14506a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 14516a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 14526a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 14536a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 14546a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 14556a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 14566a388c36SPeter Brune } 14579bd66eb0SPeter Brune } 14586a388c36SPeter Brune PetscFunctionReturn(0); 14596a388c36SPeter Brune } 14606a388c36SPeter Brune 14616f263ca3SPeter Brune /*@ 14626f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 14636f263ca3SPeter Brune 14646f263ca3SPeter Brune Input Parameters: 146578bcb3b5SPeter Brune + linesearch - linesearch context 146678bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 14676f263ca3SPeter Brune 14686f263ca3SPeter Brune Options Database Keys: 14691f8196a2SJed Brown . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch 14706f263ca3SPeter Brune 14716f263ca3SPeter Brune Notes: 14721f8196a2SJed Brown This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm. 14736f263ca3SPeter Brune 14746f263ca3SPeter Brune Level: intermediate 14756f263ca3SPeter Brune 14761a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 14776f263ca3SPeter Brune @*/ 14786f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 14796f263ca3SPeter Brune { 14806f263ca3SPeter Brune PetscFunctionBegin; 14816f263ca3SPeter Brune linesearch->norms = flg; 14826f263ca3SPeter Brune PetscFunctionReturn(0); 14836f263ca3SPeter Brune } 14846f263ca3SPeter Brune 1485f40b411bSPeter Brune /*@ 1486f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1487f40b411bSPeter Brune 1488f40b411bSPeter Brune Input Parameters: 148978bcb3b5SPeter Brune . linesearch - linesearch context 1490f40b411bSPeter Brune 1491f40b411bSPeter Brune Output Parameters: 14926232e825SPeter Brune + X - Solution vector 14936232e825SPeter Brune . F - Function vector 14946232e825SPeter Brune . Y - Search direction vector 14956232e825SPeter Brune . W - Solution work vector 14966232e825SPeter Brune - G - Function work vector 14976232e825SPeter Brune 14987bba9028SPeter Brune Notes: 14997bba9028SPeter Brune At the beginning of a line search application, X should contain a 15006232e825SPeter Brune solution and the vector F the function computed at X. At the end of the 15016232e825SPeter Brune line search application, X should contain the new solution, and F the 15026232e825SPeter Brune function evaluated at the new solution. 1503f40b411bSPeter Brune 15042a7a6963SBarry Smith These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 15052a7a6963SBarry Smith 150678bcb3b5SPeter Brune Level: advanced 1507f40b411bSPeter Brune 1508f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1509f40b411bSPeter Brune @*/ 1510bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1511bf388a1fSBarry Smith { 15126a388c36SPeter Brune PetscFunctionBegin; 1513f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15146a388c36SPeter Brune if (X) { 15156a388c36SPeter Brune PetscValidPointer(X, 2); 15166a388c36SPeter Brune *X = linesearch->vec_sol; 15176a388c36SPeter Brune } 15186a388c36SPeter Brune if (F) { 15196a388c36SPeter Brune PetscValidPointer(F, 3); 15206a388c36SPeter Brune *F = linesearch->vec_func; 15216a388c36SPeter Brune } 15226a388c36SPeter Brune if (Y) { 15236a388c36SPeter Brune PetscValidPointer(Y, 4); 15246a388c36SPeter Brune *Y = linesearch->vec_update; 15256a388c36SPeter Brune } 15266a388c36SPeter Brune if (W) { 15276a388c36SPeter Brune PetscValidPointer(W, 5); 15286a388c36SPeter Brune *W = linesearch->vec_sol_new; 15296a388c36SPeter Brune } 15306a388c36SPeter Brune if (G) { 15316a388c36SPeter Brune PetscValidPointer(G, 6); 15326a388c36SPeter Brune *G = linesearch->vec_func_new; 15336a388c36SPeter Brune } 15346a388c36SPeter Brune PetscFunctionReturn(0); 15356a388c36SPeter Brune } 15366a388c36SPeter Brune 1537f40b411bSPeter Brune /*@ 1538f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1539f40b411bSPeter Brune 1540f40b411bSPeter Brune Input Parameters: 154178bcb3b5SPeter Brune + linesearch - linesearch context 15426232e825SPeter Brune . X - Solution vector 15436232e825SPeter Brune . F - Function vector 15446232e825SPeter Brune . Y - Search direction vector 15456232e825SPeter Brune . W - Solution work vector 15466232e825SPeter Brune - G - Function work vector 1547f40b411bSPeter Brune 154878bcb3b5SPeter Brune Level: advanced 1549f40b411bSPeter Brune 1550f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1551f40b411bSPeter Brune @*/ 1552bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1553bf388a1fSBarry Smith { 15546a388c36SPeter Brune PetscFunctionBegin; 1555f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15566a388c36SPeter Brune if (X) { 15576a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 15586a388c36SPeter Brune linesearch->vec_sol = X; 15596a388c36SPeter Brune } 15606a388c36SPeter Brune if (F) { 15616a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 15626a388c36SPeter Brune linesearch->vec_func = F; 15636a388c36SPeter Brune } 15646a388c36SPeter Brune if (Y) { 15656a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 15666a388c36SPeter Brune linesearch->vec_update = Y; 15676a388c36SPeter Brune } 15686a388c36SPeter Brune if (W) { 15696a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 15706a388c36SPeter Brune linesearch->vec_sol_new = W; 15716a388c36SPeter Brune } 15726a388c36SPeter Brune if (G) { 15736a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 15746a388c36SPeter Brune linesearch->vec_func_new = G; 15756a388c36SPeter Brune } 15766a388c36SPeter Brune PetscFunctionReturn(0); 15776a388c36SPeter Brune } 15786a388c36SPeter Brune 1579e7058c64SPeter Brune /*@C 1580f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1581e7058c64SPeter Brune SNES options in the database. 1582e7058c64SPeter Brune 1583cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1584e7058c64SPeter Brune 1585e7058c64SPeter Brune Input Parameters: 1586e7058c64SPeter Brune + snes - the SNES context 1587e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1588e7058c64SPeter Brune 1589e7058c64SPeter Brune Notes: 1590e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1591e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1592e7058c64SPeter Brune 1593e7058c64SPeter Brune Level: advanced 1594e7058c64SPeter Brune 1595f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1596e7058c64SPeter Brune 1597e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1598e7058c64SPeter Brune @*/ 1599f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1600e7058c64SPeter Brune { 1601e7058c64SPeter Brune PetscErrorCode ierr; 1602e7058c64SPeter Brune 1603e7058c64SPeter Brune PetscFunctionBegin; 1604f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1605e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1606e7058c64SPeter Brune PetscFunctionReturn(0); 1607e7058c64SPeter Brune } 1608e7058c64SPeter Brune 1609e7058c64SPeter Brune /*@C 1610f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1611f1c6b773SPeter Brune SNESLineSearch options in the database. 1612e7058c64SPeter Brune 1613e7058c64SPeter Brune Not Collective 1614e7058c64SPeter Brune 1615e7058c64SPeter Brune Input Parameter: 1616cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1617e7058c64SPeter Brune 1618e7058c64SPeter Brune Output Parameter: 1619e7058c64SPeter Brune . prefix - pointer to the prefix string used 1620e7058c64SPeter Brune 16218e557f58SPeter Brune Notes: 16228e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1623e7058c64SPeter Brune sufficient length to hold the prefix. 1624e7058c64SPeter Brune 1625e7058c64SPeter Brune Level: advanced 1626e7058c64SPeter Brune 1627f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1628e7058c64SPeter Brune 1629e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1630e7058c64SPeter Brune @*/ 1631f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1632e7058c64SPeter Brune { 1633e7058c64SPeter Brune PetscErrorCode ierr; 1634e7058c64SPeter Brune 1635e7058c64SPeter Brune PetscFunctionBegin; 1636f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1637e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1638e7058c64SPeter Brune PetscFunctionReturn(0); 1639e7058c64SPeter Brune } 1640bf7f4e0aSPeter Brune 16418d359177SBarry Smith /*@C 1642fa0ddf94SBarry Smith SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1643f40b411bSPeter Brune 1644f40b411bSPeter Brune Input Parameter: 1645f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1646f40b411bSPeter Brune - nwork - the number of work vectors 1647f40b411bSPeter Brune 1648f40b411bSPeter Brune Level: developer 1649f40b411bSPeter Brune 1650fa0ddf94SBarry Smith Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1651cd7522eaSPeter Brune 1652f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1653f40b411bSPeter Brune 1654fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs() 1655f40b411bSPeter Brune @*/ 1656fa0ddf94SBarry Smith PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1657bf7f4e0aSPeter Brune { 1658bf7f4e0aSPeter Brune PetscErrorCode ierr; 1659bf388a1fSBarry Smith 1660bf7f4e0aSPeter Brune PetscFunctionBegin; 1661bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1662bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 16638d359177SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1664bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1665bf7f4e0aSPeter Brune } 1666bf7f4e0aSPeter Brune 1667f40b411bSPeter Brune /*@ 1668422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1669f40b411bSPeter Brune 1670f40b411bSPeter Brune Input Parameters: 167178bcb3b5SPeter Brune . linesearch - linesearch context 1672f40b411bSPeter Brune 1673f40b411bSPeter Brune Output Parameters: 1674422a814eSBarry Smith . result - The success or failure status 1675f40b411bSPeter Brune 1676cd7522eaSPeter Brune Notes: 1677c98378a5SBarry Smith This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1678cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1679cd7522eaSPeter Brune 1680f40b411bSPeter Brune Level: intermediate 1681f40b411bSPeter Brune 1682422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason 1683f40b411bSPeter Brune @*/ 1684422a814eSBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1685bf7f4e0aSPeter Brune { 1686bf7f4e0aSPeter Brune PetscFunctionBegin; 1687f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1688422a814eSBarry Smith PetscValidPointer(result, 2); 1689422a814eSBarry Smith *result = linesearch->result; 1690bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1691bf7f4e0aSPeter Brune } 1692bf7f4e0aSPeter Brune 1693f40b411bSPeter Brune /*@ 1694422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1695f40b411bSPeter Brune 1696f40b411bSPeter Brune Input Parameters: 169778bcb3b5SPeter Brune + linesearch - linesearch context 1698422a814eSBarry Smith - result - The success or failure status 1699f40b411bSPeter Brune 1700cd7522eaSPeter Brune Notes: 1701cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1702cd7522eaSPeter Brune the success or failure of the line search method. 1703cd7522eaSPeter Brune 170478bcb3b5SPeter Brune Level: developer 1705f40b411bSPeter Brune 1706422a814eSBarry Smith .seealso: SNESLineSearchGetSResult() 1707f40b411bSPeter Brune @*/ 1708422a814eSBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 17096a388c36SPeter Brune { 17106a388c36SPeter Brune PetscFunctionBegin; 17115fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1712422a814eSBarry Smith linesearch->result = result; 17136a388c36SPeter Brune PetscFunctionReturn(0); 17146a388c36SPeter Brune } 17156a388c36SPeter Brune 17169bd66eb0SPeter Brune /*@C 1717f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 17189bd66eb0SPeter Brune 17199bd66eb0SPeter Brune Input Parameters: 17209bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 17219bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 17229bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 17239bd66eb0SPeter Brune 17249bd66eb0SPeter Brune Logically Collective on SNES 17259bd66eb0SPeter Brune 17269bd66eb0SPeter Brune Calling sequence of projectfunc: 17279bd66eb0SPeter Brune .vb 17289bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 17299bd66eb0SPeter Brune .ve 17309bd66eb0SPeter Brune 17319bd66eb0SPeter Brune Input parameters for projectfunc: 17329bd66eb0SPeter Brune + snes - nonlinear context 17339bd66eb0SPeter Brune - X - current solution 17349bd66eb0SPeter Brune 1735cd7522eaSPeter Brune Output parameters for projectfunc: 17369bd66eb0SPeter Brune . X - Projected solution 17379bd66eb0SPeter Brune 17389bd66eb0SPeter Brune Calling sequence of normfunc: 17399bd66eb0SPeter Brune .vb 17409bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 17419bd66eb0SPeter Brune .ve 17429bd66eb0SPeter Brune 1743cd7522eaSPeter Brune Input parameters for normfunc: 17449bd66eb0SPeter Brune + snes - nonlinear context 17459bd66eb0SPeter Brune . X - current solution 17469bd66eb0SPeter Brune - F - current residual 17479bd66eb0SPeter Brune 1748cd7522eaSPeter Brune Output parameters for normfunc: 17499bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 17509bd66eb0SPeter Brune 1751cd7522eaSPeter Brune Notes: 1752cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1753cd7522eaSPeter Brune 1754cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1755cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 17569bd66eb0SPeter Brune 17579bd66eb0SPeter Brune Level: developer 17589bd66eb0SPeter Brune 17599bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 17609bd66eb0SPeter Brune 1761f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 17629bd66eb0SPeter Brune @*/ 1763f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 17649bd66eb0SPeter Brune { 17659bd66eb0SPeter Brune PetscFunctionBegin; 1766f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 17679bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 17689bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 17699bd66eb0SPeter Brune PetscFunctionReturn(0); 17709bd66eb0SPeter Brune } 17719bd66eb0SPeter Brune 17729bd66eb0SPeter Brune /*@C 1773f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 17749bd66eb0SPeter Brune 17759bd66eb0SPeter Brune Input Parameters: 1776907376e6SBarry Smith . linesearch - the line search context, obtain with SNESGetLineSearch() 17779bd66eb0SPeter Brune 17789bd66eb0SPeter Brune Output Parameters: 17799bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 17809bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 17819bd66eb0SPeter Brune 17829bd66eb0SPeter Brune Logically Collective on SNES 17839bd66eb0SPeter Brune 17849bd66eb0SPeter Brune Level: developer 17859bd66eb0SPeter Brune 17869bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 17879bd66eb0SPeter Brune 1788f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 17899bd66eb0SPeter Brune @*/ 1790f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 17919bd66eb0SPeter Brune { 17929bd66eb0SPeter Brune PetscFunctionBegin; 17939bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 17949bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 17959bd66eb0SPeter Brune PetscFunctionReturn(0); 17969bd66eb0SPeter Brune } 17979bd66eb0SPeter Brune 1798bf7f4e0aSPeter Brune /*@C 17991c84c290SBarry Smith SNESLineSearchRegister - See SNESLineSearchRegister() 1800bf7f4e0aSPeter Brune 1801bf7f4e0aSPeter Brune Level: advanced 1802bf7f4e0aSPeter Brune @*/ 1803bdf89e91SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1804bf7f4e0aSPeter Brune { 1805bf7f4e0aSPeter Brune PetscErrorCode ierr; 1806bf7f4e0aSPeter Brune 1807bf7f4e0aSPeter Brune PetscFunctionBegin; 1808a240a19fSJed Brown ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr); 1809bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1810bf7f4e0aSPeter Brune } 1811