1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2bf7f4e0aSPeter Brune 3f1c6b773SPeter Brune PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 40298fd71SBarry Smith PetscFunctionList SNESLineSearchList = NULL; 5bf7f4e0aSPeter Brune 6f1c6b773SPeter Brune PetscClassId SNESLINESEARCH_CLASSID; 7f1c6b773SPeter Brune PetscLogEvent SNESLineSearch_Apply; 8bf7f4e0aSPeter Brune 9bf7f4e0aSPeter Brune #undef __FUNCT__ 10dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorCancel" 11dcf2fd19SBarry Smith /*@ 12dcf2fd19SBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object. 13dcf2fd19SBarry Smith 14dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 15dcf2fd19SBarry Smith 16dcf2fd19SBarry Smith Input Parameters: 17dcf2fd19SBarry Smith . ls - the SNESLineSearch context 18dcf2fd19SBarry Smith 19dcf2fd19SBarry Smith Options Database Key: 20dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 21dcf2fd19SBarry Smith into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those 22dcf2fd19SBarry Smith set via the options database 23dcf2fd19SBarry Smith 24dcf2fd19SBarry Smith Notes: 25dcf2fd19SBarry Smith There is no way to clear one specific monitor from a SNESLineSearch object. 26dcf2fd19SBarry Smith 27dcf2fd19SBarry Smith This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel 28dcf2fd19SBarry Smith that one. 29dcf2fd19SBarry Smith 30dcf2fd19SBarry Smith Level: intermediate 31dcf2fd19SBarry Smith 32dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 33dcf2fd19SBarry Smith 34dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet() 35dcf2fd19SBarry Smith @*/ 36dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 37dcf2fd19SBarry Smith { 38dcf2fd19SBarry Smith PetscErrorCode ierr; 39dcf2fd19SBarry Smith PetscInt i; 40dcf2fd19SBarry Smith 41dcf2fd19SBarry Smith PetscFunctionBegin; 42dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 43dcf2fd19SBarry Smith for (i=0; i<ls->numbermonitors; i++) { 44dcf2fd19SBarry Smith if (ls->monitordestroy[i]) { 45dcf2fd19SBarry Smith ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr); 46dcf2fd19SBarry Smith } 47dcf2fd19SBarry Smith } 48dcf2fd19SBarry Smith ls->numbermonitors = 0; 49dcf2fd19SBarry Smith PetscFunctionReturn(0); 50dcf2fd19SBarry Smith } 51dcf2fd19SBarry Smith 52dcf2fd19SBarry Smith #undef __FUNCT__ 53dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitor" 54dcf2fd19SBarry Smith /*@ 55dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 56dcf2fd19SBarry Smith 57dcf2fd19SBarry Smith Collective on SNES 58dcf2fd19SBarry Smith 59dcf2fd19SBarry Smith Input Parameters: 60dcf2fd19SBarry Smith . ls - the linesearch object 61dcf2fd19SBarry Smith 62dcf2fd19SBarry Smith Notes: 63dcf2fd19SBarry Smith This routine is called by the SNES implementations. 64dcf2fd19SBarry Smith It does not typically need to be called by the user. 65dcf2fd19SBarry Smith 66dcf2fd19SBarry Smith Level: developer 67dcf2fd19SBarry Smith 68dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorSet() 69dcf2fd19SBarry Smith @*/ 70dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 71dcf2fd19SBarry Smith { 72dcf2fd19SBarry Smith PetscErrorCode ierr; 73dcf2fd19SBarry Smith PetscInt i,n = ls->numbermonitors; 74dcf2fd19SBarry Smith 75dcf2fd19SBarry Smith PetscFunctionBegin; 76dcf2fd19SBarry Smith for (i=0; i<n; i++) { 77dcf2fd19SBarry Smith ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr); 78dcf2fd19SBarry Smith } 79dcf2fd19SBarry Smith PetscFunctionReturn(0); 80dcf2fd19SBarry Smith } 81dcf2fd19SBarry Smith 82dcf2fd19SBarry Smith #undef __FUNCT__ 83dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSet" 84dcf2fd19SBarry Smith /*@C 85dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 86dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 87dcf2fd19SBarry Smith progress. 88dcf2fd19SBarry Smith 89dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 90dcf2fd19SBarry Smith 91dcf2fd19SBarry Smith Input Parameters: 92dcf2fd19SBarry Smith + ls - the SNESLineSearch context 93dcf2fd19SBarry Smith . f - the monitor function 94dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 95dcf2fd19SBarry Smith monitor routine (use NULL if no context is desired) 96dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 97dcf2fd19SBarry Smith (may be NULL) 98dcf2fd19SBarry Smith 99dcf2fd19SBarry Smith Notes: 100dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 101dcf2fd19SBarry Smith SNESLineSearchMonitorSet() multiple times; all will be called in the 102dcf2fd19SBarry Smith order in which they were set. 103dcf2fd19SBarry Smith 104dcf2fd19SBarry Smith Fortran notes: Only a single monitor function can be set for each SNESLineSearch object 105dcf2fd19SBarry Smith 106dcf2fd19SBarry Smith Level: intermediate 107dcf2fd19SBarry Smith 108dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor 109dcf2fd19SBarry Smith 110dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel() 111dcf2fd19SBarry Smith @*/ 112dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 113dcf2fd19SBarry Smith { 114dcf2fd19SBarry Smith PetscFunctionBegin; 115dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 116dcf2fd19SBarry Smith if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 117dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 118dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 119dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void*)mctx; 120dcf2fd19SBarry Smith PetscFunctionReturn(0); 121dcf2fd19SBarry Smith } 122dcf2fd19SBarry Smith 123dcf2fd19SBarry Smith #undef __FUNCT__ 124dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSolutionUpdate" 125dcf2fd19SBarry Smith /*@C 126dcf2fd19SBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries 127dcf2fd19SBarry Smith 128dcf2fd19SBarry Smith Collective on SNESLineSearch 129dcf2fd19SBarry Smith 130dcf2fd19SBarry Smith Input Parameters: 131dcf2fd19SBarry Smith + ls - the SNES linesearch object 132*d12e167eSBarry Smith - vf - the context for the monitor, in this case it is an ASCII PetscViewer and format 133dcf2fd19SBarry Smith 134dcf2fd19SBarry Smith Level: intermediate 135dcf2fd19SBarry Smith 136dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm 137dcf2fd19SBarry Smith 138dcf2fd19SBarry Smith .seealso: SNESMonitorSet(), SNESMonitorSolution() 139dcf2fd19SBarry Smith @*/ 140*d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf) 141dcf2fd19SBarry Smith { 142dcf2fd19SBarry Smith PetscErrorCode ierr; 143*d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 144dcf2fd19SBarry Smith Vec Y,W,G; 145dcf2fd19SBarry Smith 146dcf2fd19SBarry Smith PetscFunctionBegin; 147dcf2fd19SBarry Smith ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr); 148*d12e167eSBarry Smith ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 149dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr); 150dcf2fd19SBarry Smith ierr = VecView(Y,viewer);CHKERRQ(ierr); 151dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr); 152dcf2fd19SBarry Smith ierr = VecView(W,viewer);CHKERRQ(ierr); 153dcf2fd19SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr); 154dcf2fd19SBarry Smith ierr = VecView(G,viewer);CHKERRQ(ierr); 155*d12e167eSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 156dcf2fd19SBarry Smith PetscFunctionReturn(0); 157dcf2fd19SBarry Smith } 158dcf2fd19SBarry Smith 159dcf2fd19SBarry Smith #undef __FUNCT__ 160f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate" 161f40b411bSPeter Brune /*@ 162cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 163f40b411bSPeter Brune 164cd7522eaSPeter Brune Logically Collective on Comm 165f40b411bSPeter Brune 166f40b411bSPeter Brune Input Parameters: 167cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 168f40b411bSPeter Brune 169f40b411bSPeter Brune Output Parameters: 1708e557f58SPeter Brune . outlinesearch - the new linesearch context 171f40b411bSPeter Brune 172162e0bf5SPeter Brune Level: developer 173162e0bf5SPeter Brune 174162e0bf5SPeter Brune Notes: 175162e0bf5SPeter Brune The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance 176162e0bf5SPeter Brune already associated with the SNES. This function is for developer use. 177f40b411bSPeter Brune 178cd7522eaSPeter Brune .keywords: LineSearch, create, context 179f40b411bSPeter Brune 180162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch() 181f40b411bSPeter Brune @*/ 182f40b411bSPeter Brune 183bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 184bf388a1fSBarry Smith { 185bf7f4e0aSPeter Brune PetscErrorCode ierr; 186f1c6b773SPeter Brune SNESLineSearch linesearch; 187bf388a1fSBarry Smith 188bf7f4e0aSPeter Brune PetscFunctionBegin; 189ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 190f34a81feSMatthew G. Knepley ierr = SNESInitializePackage();CHKERRQ(ierr); 1910298fd71SBarry Smith *outlinesearch = NULL; 192f5af7f23SKarl Rupp 19373107ff1SLisandro Dalcin ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 194bf7f4e0aSPeter Brune 1950298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1960298fd71SBarry Smith linesearch->vec_func_new = NULL; 1970298fd71SBarry Smith linesearch->vec_sol = NULL; 1980298fd71SBarry Smith linesearch->vec_func = NULL; 1990298fd71SBarry Smith linesearch->vec_update = NULL; 2009bd66eb0SPeter Brune 201bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 202bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 203bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 204bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 205422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 206bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 207bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 208bf7f4e0aSPeter Brune linesearch->damping = 1.0; 209bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 210bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 211516fe3c3SPeter Brune linesearch->rtol = 1e-8; 212516fe3c3SPeter Brune linesearch->atol = 1e-15; 213516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2140298fd71SBarry Smith linesearch->precheckctx = NULL; 2150298fd71SBarry Smith linesearch->postcheckctx = NULL; 21659405d5eSPeter Brune linesearch->max_its = 1; 217bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 218bf7f4e0aSPeter Brune *outlinesearch = linesearch; 219bf7f4e0aSPeter Brune PetscFunctionReturn(0); 220bf7f4e0aSPeter Brune } 221bf7f4e0aSPeter Brune 222bf7f4e0aSPeter Brune #undef __FUNCT__ 223f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp" 224f40b411bSPeter Brune /*@ 22578bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 22678bcb3b5SPeter Brune any required vectors. 227f40b411bSPeter Brune 228cd7522eaSPeter Brune Collective on SNESLineSearch 229f40b411bSPeter Brune 230f40b411bSPeter Brune Input Parameters: 231f40b411bSPeter Brune . linesearch - The LineSearch instance. 232f40b411bSPeter Brune 233cd7522eaSPeter Brune Notes: 234f190f2fcSBarry Smith For most cases, this needn't be called by users or outside of SNESLineSearchApply(). 235cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 23678bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 237cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 238cd7522eaSPeter Brune allocated upfront. 239cd7522eaSPeter Brune 24078bcb3b5SPeter Brune Level: advanced 241f40b411bSPeter Brune 242f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 243f40b411bSPeter Brune 244f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 245f40b411bSPeter Brune @*/ 246f40b411bSPeter Brune 247bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 248bf388a1fSBarry Smith { 249bf7f4e0aSPeter Brune PetscErrorCode ierr; 250bf388a1fSBarry Smith 251bf7f4e0aSPeter Brune PetscFunctionBegin; 252bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 2531a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 254bf7f4e0aSPeter Brune } 255bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 2569bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 257bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 2589bd66eb0SPeter Brune } 2599bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 2609bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 2619bd66eb0SPeter Brune } 262bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 263bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 264bf7f4e0aSPeter Brune } 265ed07d7d7SPeter Brune if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);} 266bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 267bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 268bf7f4e0aSPeter Brune } 269bf7f4e0aSPeter Brune PetscFunctionReturn(0); 270bf7f4e0aSPeter Brune } 271bf7f4e0aSPeter Brune 272bf7f4e0aSPeter Brune #undef __FUNCT__ 273f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset" 274f40b411bSPeter Brune 275f40b411bSPeter Brune /*@ 276f190f2fcSBarry Smith SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search. 277f40b411bSPeter Brune 278f1c6b773SPeter Brune Collective on SNESLineSearch 279f40b411bSPeter Brune 280f40b411bSPeter Brune Input Parameters: 281f40b411bSPeter Brune . linesearch - The LineSearch instance. 282f40b411bSPeter Brune 283f190f2fcSBarry Smith Notes: Usually only called by SNESReset() 284f190f2fcSBarry Smith 285f190f2fcSBarry Smith Level: developer 286f40b411bSPeter Brune 287cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 288f40b411bSPeter Brune 289f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 290f40b411bSPeter Brune @*/ 291f40b411bSPeter Brune 292bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 293bf388a1fSBarry Smith { 294bf7f4e0aSPeter Brune PetscErrorCode ierr; 295bf388a1fSBarry Smith 296bf7f4e0aSPeter Brune PetscFunctionBegin; 297f5af7f23SKarl Rupp if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch); 298f5af7f23SKarl Rupp 299bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 300bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 301bf7f4e0aSPeter Brune 302bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 303f5af7f23SKarl Rupp 304bf7f4e0aSPeter Brune linesearch->nwork = 0; 305bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 306bf7f4e0aSPeter Brune PetscFunctionReturn(0); 307bf7f4e0aSPeter Brune } 308bf7f4e0aSPeter Brune 309ed07d7d7SPeter Brune #undef __FUNCT__ 310ed07d7d7SPeter Brune #define __FUNCT__ "SNESLineSearchSetFunction" 311ed07d7d7SPeter Brune /*@C 312f190f2fcSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search 313ed07d7d7SPeter Brune 314ed07d7d7SPeter Brune Input Parameters: 315ed07d7d7SPeter Brune . linesearch - the SNESLineSearch context 316f190f2fcSBarry Smith + func - function evaluation routine 317ed07d7d7SPeter Brune 318ed07d7d7SPeter Brune Level: developer 319ed07d7d7SPeter Brune 320f190f2fcSBarry Smith Notes: This is used internally by PETSc and not called by users 321f190f2fcSBarry Smith 322ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check 323ed07d7d7SPeter Brune 324f190f2fcSBarry Smith .seealso: SNESSetFunction() 325ed07d7d7SPeter Brune @*/ 326ed07d7d7SPeter Brune PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec)) 327ed07d7d7SPeter Brune { 328ed07d7d7SPeter Brune PetscFunctionBegin; 329ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 330ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 331ed07d7d7SPeter Brune PetscFunctionReturn(0); 332ed07d7d7SPeter Brune } 333ed07d7d7SPeter Brune 334ed07d7d7SPeter Brune 3356b2b7091SBarry Smith /*MC 336f190f2fcSBarry Smith SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called 3376b2b7091SBarry Smith 3386b2b7091SBarry Smith Synopsis: 339aaa7dc30SBarry Smith #include <petscsnes.h> 3406b2b7091SBarry Smith SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 3416b2b7091SBarry Smith 3426b2b7091SBarry Smith Input Parameters: 3436b2b7091SBarry Smith + x - solution vector 3446b2b7091SBarry Smith . y - search direction vector 3456b2b7091SBarry Smith - changed - flag to indicate the precheck changed x or y. 3466b2b7091SBarry Smith 347f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck() 348f190f2fcSBarry Smith and SNESLineSearchGetPreCheck() 349f190f2fcSBarry Smith 350878cb397SSatish Balay Level: advanced 351878cb397SSatish Balay 352f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck() 3536b2b7091SBarry Smith M*/ 3546b2b7091SBarry Smith 35586d74e61SPeter Brune #undef __FUNCT__ 356f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck" 35786d74e61SPeter Brune /*@C 358f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 359df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 360f190f2fcSBarry Smith determined the search direction. 36186d74e61SPeter Brune 362f1c6b773SPeter Brune Logically Collective on SNESLineSearch 36386d74e61SPeter Brune 36486d74e61SPeter Brune Input Parameters: 365f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 366f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence 367f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 36886d74e61SPeter Brune 36986d74e61SPeter Brune Level: intermediate 37086d74e61SPeter Brune 37186d74e61SPeter Brune .keywords: set, linesearch, pre-check 37286d74e61SPeter Brune 373f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 37486d74e61SPeter Brune @*/ 375f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx) 37686d74e61SPeter Brune { 3779bd66eb0SPeter Brune PetscFunctionBegin; 378f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 379f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 38086d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 38186d74e61SPeter Brune PetscFunctionReturn(0); 38286d74e61SPeter Brune } 38386d74e61SPeter Brune 38486d74e61SPeter Brune #undef __FUNCT__ 385f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck" 38686d74e61SPeter Brune /*@C 38753deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 38886d74e61SPeter Brune 38986d74e61SPeter Brune Input Parameters: 390f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 39186d74e61SPeter Brune 39286d74e61SPeter Brune Output Parameters: 393f190f2fcSBarry Smith + func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence 394f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 39586d74e61SPeter Brune 39686d74e61SPeter Brune Level: intermediate 39786d74e61SPeter Brune 39886d74e61SPeter Brune .keywords: get, linesearch, pre-check 39986d74e61SPeter Brune 400f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 40186d74e61SPeter Brune @*/ 402f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx) 40386d74e61SPeter Brune { 4049bd66eb0SPeter Brune PetscFunctionBegin; 405f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 406f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 40786d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 40886d74e61SPeter Brune PetscFunctionReturn(0); 40986d74e61SPeter Brune } 41086d74e61SPeter Brune 4116b2b7091SBarry Smith /*MC 412f190f2fcSBarry Smith SNESLineSearchPostCheckFunction - form of function that is called after line search is complete 4136b2b7091SBarry Smith 4146b2b7091SBarry Smith Synopsis: 415aaa7dc30SBarry Smith #include <petscsnes.h> 4166b2b7091SBarry Smith SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y, Vec w, *changed_y, PetscBool *changed_w); 4176b2b7091SBarry Smith 4186b2b7091SBarry Smith Input Parameters: 4196b2b7091SBarry Smith + x - old solution vector 4206b2b7091SBarry Smith . y - search direction vector 4216b2b7091SBarry Smith . w - new solution vector 4226b2b7091SBarry Smith . changed_y - indicates that the line search changed y 4236b2b7091SBarry Smith - changed_w - indicates that the line search changed w 4246b2b7091SBarry Smith 425f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck() 426f190f2fcSBarry Smith and SNESLineSearchGetPostCheck() 427f190f2fcSBarry Smith 428878cb397SSatish Balay Level: advanced 4296b2b7091SBarry Smith 430f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck() 4316b2b7091SBarry Smith M*/ 4326b2b7091SBarry Smith 43386d74e61SPeter Brune #undef __FUNCT__ 434f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck" 43586d74e61SPeter Brune /*@C 436f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 437f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 43886d74e61SPeter Brune 439f1c6b773SPeter Brune Logically Collective on SNESLineSearch 44086d74e61SPeter Brune 44186d74e61SPeter Brune Input Parameters: 442f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 443f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence 444f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 44586d74e61SPeter Brune 44686d74e61SPeter Brune Level: intermediate 44786d74e61SPeter Brune 44886d74e61SPeter Brune .keywords: set, linesearch, post-check 44986d74e61SPeter Brune 450f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 45186d74e61SPeter Brune @*/ 452f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 45386d74e61SPeter Brune { 45486d74e61SPeter Brune PetscFunctionBegin; 455f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 456f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 45786d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 45886d74e61SPeter Brune PetscFunctionReturn(0); 45986d74e61SPeter Brune } 46086d74e61SPeter Brune 46186d74e61SPeter Brune #undef __FUNCT__ 462f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck" 46386d74e61SPeter Brune /*@C 464f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 46586d74e61SPeter Brune 46686d74e61SPeter Brune Input Parameters: 467f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 46886d74e61SPeter Brune 46986d74e61SPeter Brune Output Parameters: 470f190f2fcSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction 471f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 47286d74e61SPeter Brune 47386d74e61SPeter Brune Level: intermediate 47486d74e61SPeter Brune 47586d74e61SPeter Brune .keywords: get, linesearch, post-check 47686d74e61SPeter Brune 477f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck() 47886d74e61SPeter Brune @*/ 479f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 48086d74e61SPeter Brune { 4819bd66eb0SPeter Brune PetscFunctionBegin; 482f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 483f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 48486d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 48586d74e61SPeter Brune PetscFunctionReturn(0); 48686d74e61SPeter Brune } 48786d74e61SPeter Brune 488bf7f4e0aSPeter Brune #undef __FUNCT__ 489f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck" 490f40b411bSPeter Brune /*@ 491f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 492f40b411bSPeter Brune 493cd7522eaSPeter Brune Logically Collective on SNESLineSearch 494f40b411bSPeter Brune 495f40b411bSPeter Brune Input Parameters: 4967b1df9c1SPeter Brune + linesearch - The linesearch instance. 4977b1df9c1SPeter Brune . X - The current solution 4987b1df9c1SPeter Brune - Y - The step direction 499f40b411bSPeter Brune 500f40b411bSPeter Brune Output Parameters: 5018e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 502f40b411bSPeter Brune 503f190f2fcSBarry Smith Level: developer 504f40b411bSPeter Brune 505f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 506f40b411bSPeter Brune 507f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 508f40b411bSPeter Brune @*/ 5097b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 510bf7f4e0aSPeter Brune { 511bf7f4e0aSPeter Brune PetscErrorCode ierr; 5125fd66863SKarl Rupp 513bf7f4e0aSPeter Brune PetscFunctionBegin; 514bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 5156b2b7091SBarry Smith if (linesearch->ops->precheck) { 5166b2b7091SBarry Smith ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr); 51738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 518bf7f4e0aSPeter Brune } 519bf7f4e0aSPeter Brune PetscFunctionReturn(0); 520bf7f4e0aSPeter Brune } 521bf7f4e0aSPeter Brune 522bf7f4e0aSPeter Brune #undef __FUNCT__ 523f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck" 524f40b411bSPeter Brune /*@ 525f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 526f40b411bSPeter Brune 527cd7522eaSPeter Brune Logically Collective on SNESLineSearch 528f40b411bSPeter Brune 529f40b411bSPeter Brune Input Parameters: 5307b1df9c1SPeter Brune + linesearch - The linesearch context 5317b1df9c1SPeter Brune . X - The last solution 5327b1df9c1SPeter Brune . Y - The step direction 5337b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 534f40b411bSPeter Brune 535f40b411bSPeter Brune Output Parameters: 53678bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 53778bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 538f40b411bSPeter Brune 539f190f2fcSBarry Smith Level: developer 540f40b411bSPeter Brune 541f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 542f40b411bSPeter Brune 543f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 544f40b411bSPeter Brune @*/ 5457b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 546bf7f4e0aSPeter Brune { 547bf7f4e0aSPeter Brune PetscErrorCode ierr; 548bf388a1fSBarry Smith 549bf7f4e0aSPeter Brune PetscFunctionBegin; 550bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 551bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 5526b2b7091SBarry Smith if (linesearch->ops->postcheck) { 5536b2b7091SBarry Smith ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr); 55438bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 55538bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 55686d74e61SPeter Brune } 55786d74e61SPeter Brune PetscFunctionReturn(0); 55886d74e61SPeter Brune } 55986d74e61SPeter Brune 56086d74e61SPeter Brune #undef __FUNCT__ 561f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard" 56286d74e61SPeter Brune /*@C 56386d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 56486d74e61SPeter Brune 565cd7522eaSPeter Brune Logically Collective on SNESLineSearch 56686d74e61SPeter Brune 56786d74e61SPeter Brune Input Arguments: 56886d74e61SPeter Brune + linesearch - linesearch context 56986d74e61SPeter Brune . X - base state for this step 57086d74e61SPeter Brune . Y - initial correction 571907376e6SBarry Smith - ctx - context for this function 57286d74e61SPeter Brune 57386d74e61SPeter Brune Output Arguments: 57486d74e61SPeter Brune + Y - correction, possibly modified 57586d74e61SPeter Brune - changed - flag indicating that Y was modified 57686d74e61SPeter Brune 57786d74e61SPeter Brune Options Database Key: 578cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 579cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 58086d74e61SPeter Brune 58186d74e61SPeter Brune Level: advanced 58286d74e61SPeter Brune 58386d74e61SPeter Brune Notes: 58486d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 58586d74e61SPeter Brune 58686d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 58786d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 58886d74e61SPeter Brune is generally not useful when using a Newton linearization. 58986d74e61SPeter Brune 59086d74e61SPeter Brune Reference: 59186d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 59286d74e61SPeter Brune 59386d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 59486d74e61SPeter Brune @*/ 595f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 59686d74e61SPeter Brune { 59786d74e61SPeter Brune PetscErrorCode ierr; 59886d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 59986d74e61SPeter Brune Vec Ylast; 60086d74e61SPeter Brune PetscScalar dot; 60186d74e61SPeter Brune PetscInt iter; 60286d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 60386d74e61SPeter Brune SNES snes; 60486d74e61SPeter Brune 60586d74e61SPeter Brune PetscFunctionBegin; 606f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 60786d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 60886d74e61SPeter Brune if (!Ylast) { 60986d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 61086d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 61186d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 61286d74e61SPeter Brune } 61386d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 61486d74e61SPeter Brune if (iter < 2) { 61586d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 61686d74e61SPeter Brune *changed = PETSC_FALSE; 61786d74e61SPeter Brune PetscFunctionReturn(0); 61886d74e61SPeter Brune } 61986d74e61SPeter Brune 62086d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 62186d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 62286d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 62386d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 624255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 62586d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 62686d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 62786d74e61SPeter Brune /* Modify the step Y */ 62886d74e61SPeter Brune PetscReal alpha,ydiffnorm; 62986d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 63086d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 63186d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 63286d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 63386d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 634c69d1a72SBarry 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); 63586d74e61SPeter Brune } else { 636c69d1a72SBarry 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); 63786d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 63886d74e61SPeter Brune *changed = PETSC_FALSE; 639bf7f4e0aSPeter Brune } 640bf7f4e0aSPeter Brune PetscFunctionReturn(0); 641bf7f4e0aSPeter Brune } 642bf7f4e0aSPeter Brune 643bf7f4e0aSPeter Brune #undef __FUNCT__ 644f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply" 645f40b411bSPeter Brune /*@ 646cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 647f40b411bSPeter Brune 648f1c6b773SPeter Brune Collective on SNESLineSearch 649f40b411bSPeter Brune 650f40b411bSPeter Brune Input Parameters: 6518e557f58SPeter Brune + linesearch - The linesearch context 6528e557f58SPeter Brune . X - The current solution 6538e557f58SPeter Brune . F - The current function 6548e557f58SPeter Brune . fnorm - The current norm 6558e557f58SPeter Brune - Y - The search direction 656f40b411bSPeter Brune 657f40b411bSPeter Brune Output Parameters: 6588e557f58SPeter Brune + X - The new solution 6598e557f58SPeter Brune . F - The new function 6608e557f58SPeter Brune - fnorm - The new function norm 661f40b411bSPeter Brune 662cd7522eaSPeter Brune Options Database Keys: 663d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell 664dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 6653c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 666cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 6673c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 6683c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 669cd7522eaSPeter Brune 670cd7522eaSPeter Brune Notes: 671cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 672cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 673cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 674cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 675cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 676cd7522eaSPeter Brune therefore may be fairly expensive. 677cd7522eaSPeter Brune 678f40b411bSPeter Brune Level: Intermediate 679f40b411bSPeter Brune 680f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 681f40b411bSPeter Brune 682cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction() 683f40b411bSPeter Brune @*/ 684bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 685bf388a1fSBarry Smith { 686bf7f4e0aSPeter Brune PetscErrorCode ierr; 687bf7f4e0aSPeter Brune 688bf388a1fSBarry Smith PetscFunctionBegin; 689f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 690bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 691bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 692bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 693bf7f4e0aSPeter Brune 694422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 695bf7f4e0aSPeter Brune 696bf7f4e0aSPeter Brune linesearch->vec_sol = X; 697bf7f4e0aSPeter Brune linesearch->vec_update = Y; 698bf7f4e0aSPeter Brune linesearch->vec_func = F; 699bf7f4e0aSPeter Brune 700f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 701bf7f4e0aSPeter Brune 702f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 703bf7f4e0aSPeter Brune 704f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 705f5af7f23SKarl Rupp else { 706bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 707bf7f4e0aSPeter Brune } 708bf7f4e0aSPeter Brune 709f1c6b773SPeter Brune ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 710bf7f4e0aSPeter Brune 711bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 712bf7f4e0aSPeter Brune 713f1c6b773SPeter Brune ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 714bf7f4e0aSPeter Brune 715f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 716bf7f4e0aSPeter Brune PetscFunctionReturn(0); 717bf7f4e0aSPeter Brune } 718bf7f4e0aSPeter Brune 719bf7f4e0aSPeter Brune #undef __FUNCT__ 720f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy" 721f40b411bSPeter Brune /*@ 722f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 723f40b411bSPeter Brune 724f1c6b773SPeter Brune Collective on SNESLineSearch 725f40b411bSPeter Brune 726f40b411bSPeter Brune Input Parameters: 7278e557f58SPeter Brune . linesearch - The linesearch context 728f40b411bSPeter Brune 729f40b411bSPeter Brune Level: Intermediate 730f40b411bSPeter Brune 73178bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 732f40b411bSPeter Brune 733cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 734f40b411bSPeter Brune @*/ 735bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 736bf388a1fSBarry Smith { 737bf7f4e0aSPeter Brune PetscErrorCode ierr; 738bf388a1fSBarry Smith 739bf7f4e0aSPeter Brune PetscFunctionBegin; 740bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 741f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 742bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 743e04113cfSBarry Smith ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr); 74422d28d08SBarry Smith ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr); 745f5af7f23SKarl Rupp if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch); 746bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 747dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr); 748e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 749bf7f4e0aSPeter Brune PetscFunctionReturn(0); 750bf7f4e0aSPeter Brune } 751bf7f4e0aSPeter Brune 752bf7f4e0aSPeter Brune #undef __FUNCT__ 753dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchSetDefaultMonitor" 754f40b411bSPeter Brune /*@ 755dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 756bf7f4e0aSPeter Brune 757bf7f4e0aSPeter Brune Input Parameters: 758dcf2fd19SBarry Smith + linesearch - the linesearch object 759dcf2fd19SBarry Smith - viewer - an ASCII PetscViewer or NULL to turn off monitor 760bf7f4e0aSPeter Brune 761dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 762bf7f4e0aSPeter Brune 763bf7f4e0aSPeter Brune Options Database: 764dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 765bf7f4e0aSPeter Brune 766bf7f4e0aSPeter Brune Level: intermediate 767bf7f4e0aSPeter Brune 768dcf2fd19SBarry Smith Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with 769*d12e167eSBarry Smith SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the 770*d12e167eSBarry Smith line search that are not visible to the other monitors. 771dcf2fd19SBarry Smith 772*d12e167eSBarry Smith .seealso: SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor() 773bf7f4e0aSPeter Brune @*/ 774dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 775bf7f4e0aSPeter Brune { 776bf7f4e0aSPeter Brune PetscErrorCode ierr; 777bf388a1fSBarry Smith 778bf7f4e0aSPeter Brune PetscFunctionBegin; 779dcf2fd19SBarry Smith if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);} 780bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 781dcf2fd19SBarry Smith linesearch->monitor = viewer; 782bf7f4e0aSPeter Brune PetscFunctionReturn(0); 783bf7f4e0aSPeter Brune } 784bf7f4e0aSPeter Brune 785bf7f4e0aSPeter Brune #undef __FUNCT__ 786dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchGetDefaultMonitor" 787f40b411bSPeter Brune /*@ 788dcf2fd19SBarry Smith SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor. 7896a388c36SPeter Brune 790f190f2fcSBarry Smith Input Parameter: 7918e557f58SPeter Brune . linesearch - linesearch context 792f40b411bSPeter Brune 793f190f2fcSBarry Smith Output Parameter: 7948e557f58SPeter Brune . monitor - monitor context 795f40b411bSPeter Brune 796f40b411bSPeter Brune Logically Collective on SNES 797f40b411bSPeter Brune 798f40b411bSPeter Brune Options Database Keys: 7998e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 800f40b411bSPeter Brune 801f40b411bSPeter Brune Level: intermediate 802f40b411bSPeter Brune 803dcf2fd19SBarry Smith .seealso: SNESLineSearchSetDefaultMonitor(), PetscViewer 804f40b411bSPeter Brune @*/ 805dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 8066a388c36SPeter Brune { 8076a388c36SPeter Brune PetscFunctionBegin; 808f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8096a388c36SPeter Brune if (monitor) { 8106a388c36SPeter Brune PetscValidPointer(monitor, 2); 8116a388c36SPeter Brune *monitor = linesearch->monitor; 8126a388c36SPeter Brune } 8136a388c36SPeter Brune PetscFunctionReturn(0); 8146a388c36SPeter Brune } 8156a388c36SPeter Brune 8166a388c36SPeter Brune #undef __FUNCT__ 817dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSetFromOptions" 818dcf2fd19SBarry Smith /*@C 819dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 820dcf2fd19SBarry Smith 821dcf2fd19SBarry Smith Collective on SNESLineSearch 822dcf2fd19SBarry Smith 823dcf2fd19SBarry Smith Input Parameters: 824dcf2fd19SBarry Smith + ls - LineSearch object you wish to monitor 825dcf2fd19SBarry Smith . name - the monitor type one is seeking 826dcf2fd19SBarry Smith . help - message indicating what monitoring is done 827dcf2fd19SBarry Smith . manual - manual page for the monitor 828dcf2fd19SBarry Smith . monitor - the monitor function 829dcf2fd19SBarry 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 830dcf2fd19SBarry Smith 831dcf2fd19SBarry Smith Level: developer 832dcf2fd19SBarry Smith 833dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 834dcf2fd19SBarry Smith PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 835dcf2fd19SBarry Smith PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 836dcf2fd19SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 837dcf2fd19SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 838dcf2fd19SBarry Smith PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 839dcf2fd19SBarry Smith PetscOptionsFList(), PetscOptionsEList() 840dcf2fd19SBarry Smith @*/ 841*d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*)) 842dcf2fd19SBarry Smith { 843dcf2fd19SBarry Smith PetscErrorCode ierr; 844dcf2fd19SBarry Smith PetscViewer viewer; 845dcf2fd19SBarry Smith PetscViewerFormat format; 846dcf2fd19SBarry Smith PetscBool flg; 847dcf2fd19SBarry Smith 848dcf2fd19SBarry Smith PetscFunctionBegin; 849dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 850dcf2fd19SBarry Smith if (flg) { 851*d12e167eSBarry Smith PetscViewerAndFormat *vf; 852*d12e167eSBarry Smith ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 853*d12e167eSBarry Smith ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 854dcf2fd19SBarry Smith if (monitorsetup) { 855*d12e167eSBarry Smith ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr); 856dcf2fd19SBarry Smith } 857*d12e167eSBarry Smith ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 858dcf2fd19SBarry Smith } 859dcf2fd19SBarry Smith PetscFunctionReturn(0); 860dcf2fd19SBarry Smith } 861dcf2fd19SBarry Smith 862dcf2fd19SBarry Smith #undef __FUNCT__ 863f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions" 864f40b411bSPeter Brune /*@ 865f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 866f40b411bSPeter Brune 867f40b411bSPeter Brune Input Parameters: 8688e557f58SPeter Brune . linesearch - linesearch context 869f40b411bSPeter Brune 870f40b411bSPeter Brune Options Database Keys: 871d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell 8723c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 8735a9b6599SPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch type 87471eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 8751a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 8761a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 8771a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 8781a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 8791a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 880dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 881dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 8828e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 883cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 8841a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 8851a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method 886f40b411bSPeter Brune 887f1c6b773SPeter Brune Logically Collective on SNESLineSearch 888f40b411bSPeter Brune 889f40b411bSPeter Brune Level: intermediate 890f40b411bSPeter Brune 8913c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard() 892f40b411bSPeter Brune @*/ 893bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 894bf388a1fSBarry Smith { 895bf7f4e0aSPeter Brune PetscErrorCode ierr; 8961a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 897bf7f4e0aSPeter Brune char type[256]; 898bf7f4e0aSPeter Brune PetscBool flg, set; 899dcf2fd19SBarry Smith PetscViewer viewer; 900bf388a1fSBarry Smith 901bf7f4e0aSPeter Brune PetscFunctionBegin; 9020f51fdf8SToby Isaac ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr); 903bf7f4e0aSPeter Brune 904bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 905f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 906a264d7a6SBarry Smith ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 907bf7f4e0aSPeter Brune if (flg) { 908f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 909bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 910f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 911bf7f4e0aSPeter Brune } 912bf7f4e0aSPeter Brune 913dcf2fd19SBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr); 914dcf2fd19SBarry Smith if (set) { 915dcf2fd19SBarry Smith ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr); 916dcf2fd19SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 917dcf2fd19SBarry Smith } 918dcf2fd19SBarry Smith ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr); 919bf7f4e0aSPeter Brune 9201a9b3a06SPeter Brune /* tolerances */ 92194ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr); 92294ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr); 92394ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr); 92494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr); 92594ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr); 92694ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr); 9277a35526eSPeter Brune 9281a9b3a06SPeter Brune /* damping parameters */ 92994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr); 9301a9b3a06SPeter Brune 93194ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr); 9321a9b3a06SPeter Brune 9331a9b3a06SPeter Brune /* precheck */ 9347a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 9357a35526eSPeter Brune if (set) { 9367a35526eSPeter Brune if (flg) { 9377a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 938f5af7f23SKarl Rupp 9397a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 9400298fd71SBarry Smith "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr); 9417a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 9427a35526eSPeter Brune } else { 9430298fd71SBarry Smith ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr); 9447a35526eSPeter Brune } 9457a35526eSPeter Brune } 94694ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr); 94794ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr); 9487a35526eSPeter Brune 9495a9b6599SPeter Brune if (linesearch->ops->setfromoptions) { 950e55864a3SBarry Smith (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr); 9515a9b6599SPeter Brune } 9525a9b6599SPeter Brune 9530633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr); 954bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 955bf7f4e0aSPeter Brune PetscFunctionReturn(0); 956bf7f4e0aSPeter Brune } 957bf7f4e0aSPeter Brune 958bf7f4e0aSPeter Brune #undef __FUNCT__ 959f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView" 960f40b411bSPeter Brune /*@ 961f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 962f40b411bSPeter Brune 963f40b411bSPeter Brune Input Parameters: 9648e557f58SPeter Brune . linesearch - linesearch context 965f40b411bSPeter Brune 966f1c6b773SPeter Brune Logically Collective on SNESLineSearch 967f40b411bSPeter Brune 968f40b411bSPeter Brune Level: intermediate 969f40b411bSPeter Brune 970dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate() 971f40b411bSPeter Brune @*/ 972bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 973bf388a1fSBarry Smith { 9747f1410a3SPeter Brune PetscErrorCode ierr; 9757f1410a3SPeter Brune PetscBool iascii; 976bf388a1fSBarry Smith 977bf7f4e0aSPeter Brune PetscFunctionBegin; 9787f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9797f1410a3SPeter Brune if (!viewer) { 980ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr); 9817f1410a3SPeter Brune } 9827f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 9837f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 984f40b411bSPeter Brune 985251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9867f1410a3SPeter Brune if (iascii) { 987dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr); 9887f1410a3SPeter Brune if (linesearch->ops->view) { 9897f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9907f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 9917f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 9927f1410a3SPeter Brune } 993c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr); 994c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr); 9957f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 9966b2b7091SBarry Smith if (linesearch->ops->precheck) { 9976b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 9987f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 9997f1410a3SPeter Brune } else { 10007f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 10017f1410a3SPeter Brune } 10027f1410a3SPeter Brune } 10036b2b7091SBarry Smith if (linesearch->ops->postcheck) { 10047f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 10057f1410a3SPeter Brune } 10067f1410a3SPeter Brune } 1007bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1008bf7f4e0aSPeter Brune } 1009bf7f4e0aSPeter Brune 1010bf7f4e0aSPeter Brune #undef __FUNCT__ 1011f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType" 1012ea5d4fccSPeter Brune /*@C 1013f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 1014f40b411bSPeter Brune 1015f190f2fcSBarry Smith Logically Collective on SNESLineSearch 1016f190f2fcSBarry Smith 1017f40b411bSPeter Brune Input Parameters: 10188e557f58SPeter Brune + linesearch - linesearch context 1019f40b411bSPeter Brune - type - The type of line search to be used 1020f40b411bSPeter Brune 1021cd7522eaSPeter Brune Available Types: 1022cd7522eaSPeter Brune + basic - Simple damping line search. 10238e557f58SPeter Brune . bt - Backtracking line search over the L2 norm of the function 10248e557f58SPeter Brune . l2 - Secant line search over the L2 norm of the function 10258e557f58SPeter Brune . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 1026d4c6564cSPatrick Farrell . nleqerr - Affine-covariant error-oriented linesearch 10278e557f58SPeter Brune - shell - User provided SNESLineSearch implementation 1028cd7522eaSPeter Brune 1029f40b411bSPeter Brune Level: intermediate 1030f40b411bSPeter Brune 1031f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 1032f40b411bSPeter Brune @*/ 103319fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 1034bf7f4e0aSPeter Brune { 1035f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 1036bf7f4e0aSPeter Brune PetscBool match; 1037bf7f4e0aSPeter Brune 1038bf7f4e0aSPeter Brune PetscFunctionBegin; 1039f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1040bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 1041bf7f4e0aSPeter Brune 1042251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 1043bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 1044bf7f4e0aSPeter Brune 10451c9cd337SJed Brown ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr); 1046bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 1047bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 1048bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 1049bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 1050f5af7f23SKarl Rupp 10510298fd71SBarry Smith linesearch->ops->destroy = NULL; 1052bf7f4e0aSPeter Brune } 1053f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 1054bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 1055bf7f4e0aSPeter Brune linesearch->ops->view = 0; 1056bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 1057bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 1058bf7f4e0aSPeter Brune 1059bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 1060bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 1061bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1062bf7f4e0aSPeter Brune } 1063bf7f4e0aSPeter Brune 1064bf7f4e0aSPeter Brune #undef __FUNCT__ 1065f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES" 1066f40b411bSPeter Brune /*@ 106778bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 1068f40b411bSPeter Brune 1069f40b411bSPeter Brune Input Parameters: 10708e557f58SPeter Brune + linesearch - linesearch context 1071f40b411bSPeter Brune - snes - The snes instance 1072f40b411bSPeter Brune 107378bcb3b5SPeter Brune Level: developer 107478bcb3b5SPeter Brune 107578bcb3b5SPeter Brune Notes: 1076f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 10777601faf0SJed Brown SNESGetLineSearch(). This routine is therefore mainly called within SNES 107878bcb3b5SPeter Brune implementations. 1079f40b411bSPeter Brune 10808141a3b9SPeter Brune Level: developer 1081f40b411bSPeter Brune 1082cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1083f40b411bSPeter Brune @*/ 1084bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1085bf388a1fSBarry Smith { 1086bf7f4e0aSPeter Brune PetscFunctionBegin; 1087f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1088bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 1089bf7f4e0aSPeter Brune linesearch->snes = snes; 1090bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1091bf7f4e0aSPeter Brune } 1092bf7f4e0aSPeter Brune 1093bf7f4e0aSPeter Brune #undef __FUNCT__ 1094f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES" 1095f40b411bSPeter Brune /*@ 10968141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 10978141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 10988141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 10998141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 1100f40b411bSPeter Brune 1101f40b411bSPeter Brune Input Parameters: 11028e557f58SPeter Brune . linesearch - linesearch context 1103f40b411bSPeter Brune 1104f40b411bSPeter Brune Output Parameters: 1105f40b411bSPeter Brune . snes - The snes instance 1106f40b411bSPeter Brune 11078141a3b9SPeter Brune Level: developer 1108f40b411bSPeter Brune 1109cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 1110f40b411bSPeter Brune @*/ 1111bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1112bf388a1fSBarry Smith { 1113bf7f4e0aSPeter Brune PetscFunctionBegin; 1114f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11156a388c36SPeter Brune PetscValidPointer(snes, 2); 1116bf7f4e0aSPeter Brune *snes = linesearch->snes; 1117bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1118bf7f4e0aSPeter Brune } 1119bf7f4e0aSPeter Brune 11206a388c36SPeter Brune #undef __FUNCT__ 1121f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda" 1122f40b411bSPeter Brune /*@ 1123f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1124f40b411bSPeter Brune 1125f40b411bSPeter Brune Input Parameters: 11268e557f58SPeter Brune . linesearch - linesearch context 1127f40b411bSPeter Brune 1128f40b411bSPeter Brune Output Parameters: 1129cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 1130f40b411bSPeter Brune 113178bcb3b5SPeter Brune Level: advanced 113278bcb3b5SPeter Brune 11338e557f58SPeter Brune Notes: 11348e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 113578bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 113678bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 113778bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 113878bcb3b5SPeter Brune 1139cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 1140f40b411bSPeter Brune @*/ 1141f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 11426a388c36SPeter Brune { 11436a388c36SPeter Brune PetscFunctionBegin; 1144f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11456a388c36SPeter Brune PetscValidPointer(lambda, 2); 11466a388c36SPeter Brune *lambda = linesearch->lambda; 11476a388c36SPeter Brune PetscFunctionReturn(0); 11486a388c36SPeter Brune } 11496a388c36SPeter Brune 11506a388c36SPeter Brune #undef __FUNCT__ 1151f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda" 1152f40b411bSPeter Brune /*@ 1153f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 1154f40b411bSPeter Brune 1155f40b411bSPeter Brune Input Parameters: 11568e557f58SPeter Brune + linesearch - linesearch context 1157f40b411bSPeter Brune - lambda - The last steplength. 1158f40b411bSPeter Brune 1159cd7522eaSPeter Brune Notes: 1160f190f2fcSBarry Smith This routine is typically used within implementations of SNESLineSearchApply() 1161cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 1162cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1163cd7522eaSPeter Brune as an inner scaling parameter. 1164cd7522eaSPeter Brune 116578bcb3b5SPeter Brune Level: advanced 1166f40b411bSPeter Brune 1167f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 1168f40b411bSPeter Brune @*/ 1169f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 11706a388c36SPeter Brune { 11716a388c36SPeter Brune PetscFunctionBegin; 1172f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11736a388c36SPeter Brune linesearch->lambda = lambda; 11746a388c36SPeter Brune PetscFunctionReturn(0); 11756a388c36SPeter Brune } 11766a388c36SPeter Brune 11776a388c36SPeter Brune #undef __FUNCT__ 1178f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances" 1179f40b411bSPeter Brune /*@ 11803c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 118178bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 118278bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 118378bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1184f40b411bSPeter Brune 1185f40b411bSPeter Brune Input Parameters: 11868e557f58SPeter Brune . linesearch - linesearch context 1187f40b411bSPeter Brune 1188f40b411bSPeter Brune Output Parameters: 1189516fe3c3SPeter Brune + steptol - The minimum steplength 11906cc8e53bSPeter Brune . maxstep - The maximum steplength 1191516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1192516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1193516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1194516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1195f40b411bSPeter Brune 119678bcb3b5SPeter Brune Level: intermediate 119778bcb3b5SPeter Brune 119878bcb3b5SPeter Brune Notes: 119978bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 12003c7d6663SPeter Brune the type requires. 1201516fe3c3SPeter Brune 1202f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 1203f40b411bSPeter Brune @*/ 1204f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 12056a388c36SPeter Brune { 12066a388c36SPeter Brune PetscFunctionBegin; 1207f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1208516fe3c3SPeter Brune if (steptol) { 12096a388c36SPeter Brune PetscValidPointer(steptol, 2); 12106a388c36SPeter Brune *steptol = linesearch->steptol; 1211516fe3c3SPeter Brune } 1212516fe3c3SPeter Brune if (maxstep) { 1213516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 1214516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1215516fe3c3SPeter Brune } 1216516fe3c3SPeter Brune if (rtol) { 1217516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 1218516fe3c3SPeter Brune *rtol = linesearch->rtol; 1219516fe3c3SPeter Brune } 1220516fe3c3SPeter Brune if (atol) { 1221516fe3c3SPeter Brune PetscValidPointer(atol, 5); 1222516fe3c3SPeter Brune *atol = linesearch->atol; 1223516fe3c3SPeter Brune } 1224516fe3c3SPeter Brune if (ltol) { 1225516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 1226516fe3c3SPeter Brune *ltol = linesearch->ltol; 1227516fe3c3SPeter Brune } 1228516fe3c3SPeter Brune if (max_its) { 1229516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 1230516fe3c3SPeter Brune *max_its = linesearch->max_its; 1231516fe3c3SPeter Brune } 12326a388c36SPeter Brune PetscFunctionReturn(0); 12336a388c36SPeter Brune } 12346a388c36SPeter Brune 12356a388c36SPeter Brune #undef __FUNCT__ 1236f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances" 1237f40b411bSPeter Brune /*@ 12383c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 123978bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 124078bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 124178bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1242f40b411bSPeter Brune 1243f40b411bSPeter Brune Input Parameters: 12448e557f58SPeter Brune + linesearch - linesearch context 1245516fe3c3SPeter Brune . steptol - The minimum steplength 12466cc8e53bSPeter Brune . maxstep - The maximum steplength 1247516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1248516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1249516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1250516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1251f40b411bSPeter Brune 125278bcb3b5SPeter Brune Notes: 12533c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 125478bcb3b5SPeter Brune place of an argument. 1255f40b411bSPeter Brune 125678bcb3b5SPeter Brune Level: intermediate 1257516fe3c3SPeter Brune 1258f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 1259f40b411bSPeter Brune @*/ 1260f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 12616a388c36SPeter Brune { 12626a388c36SPeter Brune PetscFunctionBegin; 1263f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1264d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1265d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1266d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1267d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1268d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1269d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1270d3952184SSatish Balay 1271d3952184SSatish Balay if (steptol!= PETSC_DEFAULT) { 1272ce94432eSBarry Smith if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 12736a388c36SPeter Brune linesearch->steptol = steptol; 1274d3952184SSatish Balay } 1275d3952184SSatish Balay 1276d3952184SSatish Balay if (maxstep!= PETSC_DEFAULT) { 1277ce94432eSBarry Smith if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1278516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1279d3952184SSatish Balay } 1280d3952184SSatish Balay 1281d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1282ce94432eSBarry 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); 1283516fe3c3SPeter Brune linesearch->rtol = rtol; 1284d3952184SSatish Balay } 1285d3952184SSatish Balay 1286d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1287ce94432eSBarry Smith if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1288516fe3c3SPeter Brune linesearch->atol = atol; 1289d3952184SSatish Balay } 1290d3952184SSatish Balay 1291d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1292ce94432eSBarry Smith if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1293516fe3c3SPeter Brune linesearch->ltol = ltol; 1294d3952184SSatish Balay } 1295d3952184SSatish Balay 1296d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1297ce94432eSBarry Smith if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1298516fe3c3SPeter Brune linesearch->max_its = max_its; 1299d3952184SSatish Balay } 13006a388c36SPeter Brune PetscFunctionReturn(0); 13016a388c36SPeter Brune } 13026a388c36SPeter Brune 13036a388c36SPeter Brune #undef __FUNCT__ 1304f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping" 1305f40b411bSPeter Brune /*@ 1306f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1307f40b411bSPeter Brune 1308f40b411bSPeter Brune Input Parameters: 13098e557f58SPeter Brune . linesearch - linesearch context 1310f40b411bSPeter Brune 1311f40b411bSPeter Brune Output Parameters: 13128e557f58SPeter Brune . damping - The damping parameter 1313f40b411bSPeter Brune 131478bcb3b5SPeter Brune Level: advanced 1315f40b411bSPeter Brune 131678bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1317f40b411bSPeter Brune @*/ 1318f40b411bSPeter Brune 1319f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 13206a388c36SPeter Brune { 13216a388c36SPeter Brune PetscFunctionBegin; 1322f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13236a388c36SPeter Brune PetscValidPointer(damping, 2); 13246a388c36SPeter Brune *damping = linesearch->damping; 13256a388c36SPeter Brune PetscFunctionReturn(0); 13266a388c36SPeter Brune } 13276a388c36SPeter Brune 13286a388c36SPeter Brune #undef __FUNCT__ 1329f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping" 1330f40b411bSPeter Brune /*@ 1331f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1332f40b411bSPeter Brune 1333f40b411bSPeter Brune Input Parameters: 133403fd4120SBarry Smith + linesearch - linesearch context 133503fd4120SBarry Smith - damping - The damping parameter 1336f40b411bSPeter Brune 133703fd4120SBarry Smith Options Database: 133803fd4120SBarry Smith . -snes_linesearch_damping 1339f40b411bSPeter Brune Level: intermediate 1340f40b411bSPeter Brune 1341cd7522eaSPeter Brune Notes: 1342cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1343cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 134478bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1345cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1346cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1347cd7522eaSPeter Brune 1348f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1349f40b411bSPeter Brune @*/ 1350f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 13516a388c36SPeter Brune { 13526a388c36SPeter Brune PetscFunctionBegin; 1353f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13546a388c36SPeter Brune linesearch->damping = damping; 13556a388c36SPeter Brune PetscFunctionReturn(0); 13566a388c36SPeter Brune } 13576a388c36SPeter Brune 13586a388c36SPeter Brune #undef __FUNCT__ 135959405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder" 136059405d5eSPeter Brune /*@ 136159405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 136259405d5eSPeter Brune 136359405d5eSPeter Brune Input Parameters: 136478bcb3b5SPeter Brune . linesearch - linesearch context 136559405d5eSPeter Brune 136659405d5eSPeter Brune Output Parameters: 13678e557f58SPeter Brune . order - The order 136859405d5eSPeter Brune 136978bcb3b5SPeter Brune Possible Values for order: 13703c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 13713c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 13723c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 137378bcb3b5SPeter Brune 137459405d5eSPeter Brune Level: intermediate 137559405d5eSPeter Brune 137659405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 137759405d5eSPeter Brune @*/ 137859405d5eSPeter Brune 1379b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 138059405d5eSPeter Brune { 138159405d5eSPeter Brune PetscFunctionBegin; 138259405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 138359405d5eSPeter Brune PetscValidPointer(order, 2); 138459405d5eSPeter Brune *order = linesearch->order; 138559405d5eSPeter Brune PetscFunctionReturn(0); 138659405d5eSPeter Brune } 138759405d5eSPeter Brune 138859405d5eSPeter Brune #undef __FUNCT__ 138959405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder" 139059405d5eSPeter Brune /*@ 139159405d5eSPeter Brune SNESLineSearchSetOrder - Sets the line search damping paramter. 139259405d5eSPeter Brune 139359405d5eSPeter Brune Input Parameters: 139478bcb3b5SPeter Brune . linesearch - linesearch context 139578bcb3b5SPeter Brune . order - The damping parameter 139659405d5eSPeter Brune 139759405d5eSPeter Brune Level: intermediate 139859405d5eSPeter Brune 139978bcb3b5SPeter Brune Possible Values for order: 14003c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 14013c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 14023c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 140378bcb3b5SPeter Brune 140459405d5eSPeter Brune Notes: 140559405d5eSPeter Brune Variable orders are supported by the following line searches: 140678bcb3b5SPeter Brune + bt - cubic and quadratic 140778bcb3b5SPeter Brune - cp - linear and quadratic 140859405d5eSPeter Brune 140959405d5eSPeter Brune .seealso: SNESLineSearchGetOrder() 141059405d5eSPeter Brune @*/ 1411b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 141259405d5eSPeter Brune { 141359405d5eSPeter Brune PetscFunctionBegin; 141459405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 141559405d5eSPeter Brune linesearch->order = order; 141659405d5eSPeter Brune PetscFunctionReturn(0); 141759405d5eSPeter Brune } 141859405d5eSPeter Brune 141959405d5eSPeter Brune #undef __FUNCT__ 1420f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms" 1421f40b411bSPeter Brune /*@ 1422f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1423f40b411bSPeter Brune 1424f40b411bSPeter Brune Input Parameters: 142578bcb3b5SPeter Brune . linesearch - linesearch context 1426f40b411bSPeter Brune 1427f40b411bSPeter Brune Output Parameters: 1428f40b411bSPeter Brune + xnorm - The norm of the current solution 1429f40b411bSPeter Brune . fnorm - The norm of the current function 1430f40b411bSPeter Brune - ynorm - The norm of the current update 1431f40b411bSPeter Brune 1432cd7522eaSPeter Brune Notes: 1433cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1434cd7522eaSPeter Brune 143578bcb3b5SPeter Brune Level: developer 1436f40b411bSPeter Brune 1437f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1438f40b411bSPeter Brune @*/ 1439f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1440bf7f4e0aSPeter Brune { 1441bf7f4e0aSPeter Brune PetscFunctionBegin; 1442f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1443f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1444f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1445f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 1446bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1447bf7f4e0aSPeter Brune } 1448bf7f4e0aSPeter Brune 1449e7058c64SPeter Brune #undef __FUNCT__ 1450f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms" 1451f40b411bSPeter Brune /*@ 1452f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1453f40b411bSPeter Brune 1454f40b411bSPeter Brune Input Parameters: 145578bcb3b5SPeter Brune + linesearch - linesearch context 1456f40b411bSPeter Brune . xnorm - The norm of the current solution 1457f40b411bSPeter Brune . fnorm - The norm of the current function 1458f40b411bSPeter Brune - ynorm - The norm of the current update 1459f40b411bSPeter Brune 146078bcb3b5SPeter Brune Level: advanced 1461f40b411bSPeter Brune 1462f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1463f40b411bSPeter Brune @*/ 1464f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 14656a388c36SPeter Brune { 14666a388c36SPeter Brune PetscFunctionBegin; 1467f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14686a388c36SPeter Brune linesearch->xnorm = xnorm; 14696a388c36SPeter Brune linesearch->fnorm = fnorm; 14706a388c36SPeter Brune linesearch->ynorm = ynorm; 14716a388c36SPeter Brune PetscFunctionReturn(0); 14726a388c36SPeter Brune } 14736a388c36SPeter Brune 14746a388c36SPeter Brune #undef __FUNCT__ 1475f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms" 1476f40b411bSPeter Brune /*@ 1477f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1478f40b411bSPeter Brune 1479f40b411bSPeter Brune Input Parameters: 148078bcb3b5SPeter Brune . linesearch - linesearch context 1481f40b411bSPeter Brune 1482f40b411bSPeter Brune Options Database Keys: 14838e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1484f40b411bSPeter Brune 1485f40b411bSPeter Brune Level: intermediate 1486f40b411bSPeter Brune 148778bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1488f40b411bSPeter Brune @*/ 1489f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 14906a388c36SPeter Brune { 14916a388c36SPeter Brune PetscErrorCode ierr; 14929bd66eb0SPeter Brune SNES snes; 1493bf388a1fSBarry Smith 14946a388c36SPeter Brune PetscFunctionBegin; 14956a388c36SPeter Brune if (linesearch->norms) { 14969bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1497f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 14989bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 14999bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15009bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 15019bd66eb0SPeter Brune } else { 15026a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 15036a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 15046a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15056a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 15066a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 15076a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 15086a388c36SPeter Brune } 15099bd66eb0SPeter Brune } 15106a388c36SPeter Brune PetscFunctionReturn(0); 15116a388c36SPeter Brune } 15126a388c36SPeter Brune 15136f263ca3SPeter Brune #undef __FUNCT__ 15146f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms" 15156f263ca3SPeter Brune /*@ 15166f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 15176f263ca3SPeter Brune 15186f263ca3SPeter Brune Input Parameters: 151978bcb3b5SPeter Brune + linesearch - linesearch context 152078bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 15216f263ca3SPeter Brune 15226f263ca3SPeter Brune Options Database Keys: 15238e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 15246f263ca3SPeter Brune 15256f263ca3SPeter Brune Notes: 15261a4f838cSPeter Brune This is most relevant to the SNESLINESEARCHBASIC line search type. 15276f263ca3SPeter Brune 15286f263ca3SPeter Brune Level: intermediate 15296f263ca3SPeter Brune 15301a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 15316f263ca3SPeter Brune @*/ 15326f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 15336f263ca3SPeter Brune { 15346f263ca3SPeter Brune PetscFunctionBegin; 15356f263ca3SPeter Brune linesearch->norms = flg; 15366f263ca3SPeter Brune PetscFunctionReturn(0); 15376f263ca3SPeter Brune } 15386f263ca3SPeter Brune 15396a388c36SPeter Brune #undef __FUNCT__ 1540f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs" 1541f40b411bSPeter Brune /*@ 1542f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1543f40b411bSPeter Brune 1544f40b411bSPeter Brune Input Parameters: 154578bcb3b5SPeter Brune . linesearch - linesearch context 1546f40b411bSPeter Brune 1547f40b411bSPeter Brune Output Parameters: 15486232e825SPeter Brune + X - Solution vector 15496232e825SPeter Brune . F - Function vector 15506232e825SPeter Brune . Y - Search direction vector 15516232e825SPeter Brune . W - Solution work vector 15526232e825SPeter Brune - G - Function work vector 15536232e825SPeter Brune 15547bba9028SPeter Brune Notes: 15557bba9028SPeter Brune At the beginning of a line search application, X should contain a 15566232e825SPeter Brune solution and the vector F the function computed at X. At the end of the 15576232e825SPeter Brune line search application, X should contain the new solution, and F the 15586232e825SPeter Brune function evaluated at the new solution. 1559f40b411bSPeter Brune 15602a7a6963SBarry Smith These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 15612a7a6963SBarry Smith 156278bcb3b5SPeter Brune Level: advanced 1563f40b411bSPeter Brune 1564f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1565f40b411bSPeter Brune @*/ 1566bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1567bf388a1fSBarry Smith { 15686a388c36SPeter Brune PetscFunctionBegin; 1569f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15706a388c36SPeter Brune if (X) { 15716a388c36SPeter Brune PetscValidPointer(X, 2); 15726a388c36SPeter Brune *X = linesearch->vec_sol; 15736a388c36SPeter Brune } 15746a388c36SPeter Brune if (F) { 15756a388c36SPeter Brune PetscValidPointer(F, 3); 15766a388c36SPeter Brune *F = linesearch->vec_func; 15776a388c36SPeter Brune } 15786a388c36SPeter Brune if (Y) { 15796a388c36SPeter Brune PetscValidPointer(Y, 4); 15806a388c36SPeter Brune *Y = linesearch->vec_update; 15816a388c36SPeter Brune } 15826a388c36SPeter Brune if (W) { 15836a388c36SPeter Brune PetscValidPointer(W, 5); 15846a388c36SPeter Brune *W = linesearch->vec_sol_new; 15856a388c36SPeter Brune } 15866a388c36SPeter Brune if (G) { 15876a388c36SPeter Brune PetscValidPointer(G, 6); 15886a388c36SPeter Brune *G = linesearch->vec_func_new; 15896a388c36SPeter Brune } 15906a388c36SPeter Brune PetscFunctionReturn(0); 15916a388c36SPeter Brune } 15926a388c36SPeter Brune 15936a388c36SPeter Brune #undef __FUNCT__ 1594f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs" 1595f40b411bSPeter Brune /*@ 1596f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1597f40b411bSPeter Brune 1598f40b411bSPeter Brune Input Parameters: 159978bcb3b5SPeter Brune + linesearch - linesearch context 16006232e825SPeter Brune . X - Solution vector 16016232e825SPeter Brune . F - Function vector 16026232e825SPeter Brune . Y - Search direction vector 16036232e825SPeter Brune . W - Solution work vector 16046232e825SPeter Brune - G - Function work vector 1605f40b411bSPeter Brune 160678bcb3b5SPeter Brune Level: advanced 1607f40b411bSPeter Brune 1608f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1609f40b411bSPeter Brune @*/ 1610bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1611bf388a1fSBarry Smith { 16126a388c36SPeter Brune PetscFunctionBegin; 1613f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 16146a388c36SPeter Brune if (X) { 16156a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 16166a388c36SPeter Brune linesearch->vec_sol = X; 16176a388c36SPeter Brune } 16186a388c36SPeter Brune if (F) { 16196a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 16206a388c36SPeter Brune linesearch->vec_func = F; 16216a388c36SPeter Brune } 16226a388c36SPeter Brune if (Y) { 16236a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 16246a388c36SPeter Brune linesearch->vec_update = Y; 16256a388c36SPeter Brune } 16266a388c36SPeter Brune if (W) { 16276a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 16286a388c36SPeter Brune linesearch->vec_sol_new = W; 16296a388c36SPeter Brune } 16306a388c36SPeter Brune if (G) { 16316a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 16326a388c36SPeter Brune linesearch->vec_func_new = G; 16336a388c36SPeter Brune } 16346a388c36SPeter Brune PetscFunctionReturn(0); 16356a388c36SPeter Brune } 16366a388c36SPeter Brune 16376a388c36SPeter Brune #undef __FUNCT__ 1638f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1639e7058c64SPeter Brune /*@C 1640f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1641e7058c64SPeter Brune SNES options in the database. 1642e7058c64SPeter Brune 1643cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1644e7058c64SPeter Brune 1645e7058c64SPeter Brune Input Parameters: 1646e7058c64SPeter Brune + snes - the SNES context 1647e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1648e7058c64SPeter Brune 1649e7058c64SPeter Brune Notes: 1650e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1651e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1652e7058c64SPeter Brune 1653e7058c64SPeter Brune Level: advanced 1654e7058c64SPeter Brune 1655f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1656e7058c64SPeter Brune 1657e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1658e7058c64SPeter Brune @*/ 1659f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1660e7058c64SPeter Brune { 1661e7058c64SPeter Brune PetscErrorCode ierr; 1662e7058c64SPeter Brune 1663e7058c64SPeter Brune PetscFunctionBegin; 1664f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1665e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1666e7058c64SPeter Brune PetscFunctionReturn(0); 1667e7058c64SPeter Brune } 1668e7058c64SPeter Brune 1669e7058c64SPeter Brune #undef __FUNCT__ 1670f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1671e7058c64SPeter Brune /*@C 1672f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1673f1c6b773SPeter Brune SNESLineSearch options in the database. 1674e7058c64SPeter Brune 1675e7058c64SPeter Brune Not Collective 1676e7058c64SPeter Brune 1677e7058c64SPeter Brune Input Parameter: 1678cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1679e7058c64SPeter Brune 1680e7058c64SPeter Brune Output Parameter: 1681e7058c64SPeter Brune . prefix - pointer to the prefix string used 1682e7058c64SPeter Brune 16838e557f58SPeter Brune Notes: 16848e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1685e7058c64SPeter Brune sufficient length to hold the prefix. 1686e7058c64SPeter Brune 1687e7058c64SPeter Brune Level: advanced 1688e7058c64SPeter Brune 1689f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1690e7058c64SPeter Brune 1691e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1692e7058c64SPeter Brune @*/ 1693f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1694e7058c64SPeter Brune { 1695e7058c64SPeter Brune PetscErrorCode ierr; 1696e7058c64SPeter Brune 1697e7058c64SPeter Brune PetscFunctionBegin; 1698f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1699e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1700e7058c64SPeter Brune PetscFunctionReturn(0); 1701e7058c64SPeter Brune } 1702bf7f4e0aSPeter Brune 1703bf7f4e0aSPeter Brune #undef __FUNCT__ 1704fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs" 17058d359177SBarry Smith /*@C 1706fa0ddf94SBarry Smith SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1707f40b411bSPeter Brune 1708f40b411bSPeter Brune Input Parameter: 1709f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1710f40b411bSPeter Brune - nwork - the number of work vectors 1711f40b411bSPeter Brune 1712f40b411bSPeter Brune Level: developer 1713f40b411bSPeter Brune 1714fa0ddf94SBarry Smith Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1715cd7522eaSPeter Brune 1716f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1717f40b411bSPeter Brune 1718fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs() 1719f40b411bSPeter Brune @*/ 1720fa0ddf94SBarry Smith PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1721bf7f4e0aSPeter Brune { 1722bf7f4e0aSPeter Brune PetscErrorCode ierr; 1723bf388a1fSBarry Smith 1724bf7f4e0aSPeter Brune PetscFunctionBegin; 1725bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1726bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 17278d359177SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1728bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1729bf7f4e0aSPeter Brune } 1730bf7f4e0aSPeter Brune 1731bf7f4e0aSPeter Brune #undef __FUNCT__ 1732422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchGetReason" 1733f40b411bSPeter Brune /*@ 1734422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1735f40b411bSPeter Brune 1736f40b411bSPeter Brune Input Parameters: 173778bcb3b5SPeter Brune . linesearch - linesearch context 1738f40b411bSPeter Brune 1739f40b411bSPeter Brune Output Parameters: 1740422a814eSBarry Smith . result - The success or failure status 1741f40b411bSPeter Brune 1742cd7522eaSPeter Brune Notes: 1743c98378a5SBarry Smith This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1744cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1745cd7522eaSPeter Brune 1746f40b411bSPeter Brune Level: intermediate 1747f40b411bSPeter Brune 1748422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason 1749f40b411bSPeter Brune @*/ 1750422a814eSBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1751bf7f4e0aSPeter Brune { 1752bf7f4e0aSPeter Brune PetscFunctionBegin; 1753f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1754422a814eSBarry Smith PetscValidPointer(result, 2); 1755422a814eSBarry Smith *result = linesearch->result; 1756bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1757bf7f4e0aSPeter Brune } 1758bf7f4e0aSPeter Brune 1759bf7f4e0aSPeter Brune #undef __FUNCT__ 1760422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchSetReason" 1761f40b411bSPeter Brune /*@ 1762422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1763f40b411bSPeter Brune 1764f40b411bSPeter Brune Input Parameters: 176578bcb3b5SPeter Brune + linesearch - linesearch context 1766422a814eSBarry Smith - result - The success or failure status 1767f40b411bSPeter Brune 1768cd7522eaSPeter Brune Notes: 1769cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1770cd7522eaSPeter Brune the success or failure of the line search method. 1771cd7522eaSPeter Brune 177278bcb3b5SPeter Brune Level: developer 1773f40b411bSPeter Brune 1774422a814eSBarry Smith .seealso: SNESLineSearchGetSResult() 1775f40b411bSPeter Brune @*/ 1776422a814eSBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 17776a388c36SPeter Brune { 17786a388c36SPeter Brune PetscFunctionBegin; 17795fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1780422a814eSBarry Smith linesearch->result = result; 17816a388c36SPeter Brune PetscFunctionReturn(0); 17826a388c36SPeter Brune } 17836a388c36SPeter Brune 17846a388c36SPeter Brune #undef __FUNCT__ 1785f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions" 17869bd66eb0SPeter Brune /*@C 1787f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 17889bd66eb0SPeter Brune 17899bd66eb0SPeter Brune Input Parameters: 17909bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 17919bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 17929bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 17939bd66eb0SPeter Brune 17949bd66eb0SPeter Brune Logically Collective on SNES 17959bd66eb0SPeter Brune 17969bd66eb0SPeter Brune Calling sequence of projectfunc: 17979bd66eb0SPeter Brune .vb 17989bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 17999bd66eb0SPeter Brune .ve 18009bd66eb0SPeter Brune 18019bd66eb0SPeter Brune Input parameters for projectfunc: 18029bd66eb0SPeter Brune + snes - nonlinear context 18039bd66eb0SPeter Brune - X - current solution 18049bd66eb0SPeter Brune 1805cd7522eaSPeter Brune Output parameters for projectfunc: 18069bd66eb0SPeter Brune . X - Projected solution 18079bd66eb0SPeter Brune 18089bd66eb0SPeter Brune Calling sequence of normfunc: 18099bd66eb0SPeter Brune .vb 18109bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 18119bd66eb0SPeter Brune .ve 18129bd66eb0SPeter Brune 1813cd7522eaSPeter Brune Input parameters for normfunc: 18149bd66eb0SPeter Brune + snes - nonlinear context 18159bd66eb0SPeter Brune . X - current solution 18169bd66eb0SPeter Brune - F - current residual 18179bd66eb0SPeter Brune 1818cd7522eaSPeter Brune Output parameters for normfunc: 18199bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 18209bd66eb0SPeter Brune 1821cd7522eaSPeter Brune Notes: 1822cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1823cd7522eaSPeter Brune 1824cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1825cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 18269bd66eb0SPeter Brune 18279bd66eb0SPeter Brune Level: developer 18289bd66eb0SPeter Brune 18299bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 18309bd66eb0SPeter Brune 1831f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 18329bd66eb0SPeter Brune @*/ 1833f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 18349bd66eb0SPeter Brune { 18359bd66eb0SPeter Brune PetscFunctionBegin; 1836f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 18379bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 18389bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 18399bd66eb0SPeter Brune PetscFunctionReturn(0); 18409bd66eb0SPeter Brune } 18419bd66eb0SPeter Brune 18429bd66eb0SPeter Brune #undef __FUNCT__ 1843f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions" 18449bd66eb0SPeter Brune /*@C 1845f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 18469bd66eb0SPeter Brune 18479bd66eb0SPeter Brune Input Parameters: 1848907376e6SBarry Smith . linesearch - the line search context, obtain with SNESGetLineSearch() 18499bd66eb0SPeter Brune 18509bd66eb0SPeter Brune Output Parameters: 18519bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 18529bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 18539bd66eb0SPeter Brune 18549bd66eb0SPeter Brune Logically Collective on SNES 18559bd66eb0SPeter Brune 18569bd66eb0SPeter Brune Level: developer 18579bd66eb0SPeter Brune 18589bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 18599bd66eb0SPeter Brune 1860f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 18619bd66eb0SPeter Brune @*/ 1862f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 18639bd66eb0SPeter Brune { 18649bd66eb0SPeter Brune PetscFunctionBegin; 18659bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 18669bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 18679bd66eb0SPeter Brune PetscFunctionReturn(0); 18689bd66eb0SPeter Brune } 18699bd66eb0SPeter Brune 18709bd66eb0SPeter Brune #undef __FUNCT__ 1871f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister" 1872bf7f4e0aSPeter Brune /*@C 18731c84c290SBarry Smith SNESLineSearchRegister - See SNESLineSearchRegister() 1874bf7f4e0aSPeter Brune 1875bf7f4e0aSPeter Brune Level: advanced 1876bf7f4e0aSPeter Brune @*/ 1877bdf89e91SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1878bf7f4e0aSPeter Brune { 1879bf7f4e0aSPeter Brune PetscErrorCode ierr; 1880bf7f4e0aSPeter Brune 1881bf7f4e0aSPeter Brune PetscFunctionBegin; 1882a240a19fSJed Brown ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr); 1883bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1884bf7f4e0aSPeter Brune } 1885