1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2bf7f4e0aSPeter Brune 3f1c6b773SPeter Brune PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 40298fd71SBarry Smith PetscFunctionList SNESLineSearchList = NULL; 5bf7f4e0aSPeter Brune 6f1c6b773SPeter Brune PetscClassId SNESLINESEARCH_CLASSID; 757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply; 8bf7f4e0aSPeter Brune 9dcf2fd19SBarry Smith /*@ 10f6dfbefdSBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object. 11dcf2fd19SBarry Smith 12c3339decSBarry Smith Logically Collective 13dcf2fd19SBarry Smith 14dcf2fd19SBarry Smith Input Parameters: 15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context 16dcf2fd19SBarry Smith 17dcf2fd19SBarry Smith Options Database Key: 18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 19f6dfbefdSBarry Smith into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those 20dcf2fd19SBarry Smith set via the options database 21dcf2fd19SBarry Smith 22*20f4b53cSBarry Smith Level: advanced 23*20f4b53cSBarry Smith 24dcf2fd19SBarry Smith Notes: 25f6dfbefdSBarry Smith There is no way to clear one specific monitor from a `SNESLineSearch` object. 26dcf2fd19SBarry Smith 27f6dfbefdSBarry Smith This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel 28dcf2fd19SBarry Smith that one. 29dcf2fd19SBarry Smith 30f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()` 31dcf2fd19SBarry Smith @*/ 32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 33d71ae5a4SJacob Faibussowitsch { 34dcf2fd19SBarry Smith PetscInt i; 35dcf2fd19SBarry Smith 36dcf2fd19SBarry Smith PetscFunctionBegin; 37dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 38dcf2fd19SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 3948a46eb9SPierre Jolivet if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i])); 40dcf2fd19SBarry Smith } 41dcf2fd19SBarry Smith ls->numbermonitors = 0; 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43dcf2fd19SBarry Smith } 44dcf2fd19SBarry Smith 45dcf2fd19SBarry Smith /*@ 46dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 47dcf2fd19SBarry Smith 48c3339decSBarry Smith Collective 49dcf2fd19SBarry Smith 50dcf2fd19SBarry Smith Input Parameters: 51dcf2fd19SBarry Smith . ls - the linesearch object 52dcf2fd19SBarry Smith 53f6dfbefdSBarry Smith Note: 54*20f4b53cSBarry Smith This routine is called by the `SNES` implementations. 55dcf2fd19SBarry Smith It does not typically need to be called by the user. 56dcf2fd19SBarry Smith 57dcf2fd19SBarry Smith Level: developer 58dcf2fd19SBarry Smith 59f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()` 60dcf2fd19SBarry Smith @*/ 61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 62d71ae5a4SJacob Faibussowitsch { 63dcf2fd19SBarry Smith PetscInt i, n = ls->numbermonitors; 64dcf2fd19SBarry Smith 65dcf2fd19SBarry Smith PetscFunctionBegin; 6648a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i])); 673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68dcf2fd19SBarry Smith } 69dcf2fd19SBarry Smith 70dcf2fd19SBarry Smith /*@C 71dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 72dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 73dcf2fd19SBarry Smith progress. 74dcf2fd19SBarry Smith 75c3339decSBarry Smith Logically Collective 76dcf2fd19SBarry Smith 77dcf2fd19SBarry Smith Input Parameters: 78f6dfbefdSBarry Smith + ls - the `SNESLineSearch` context 79dcf2fd19SBarry Smith . f - the monitor function 80dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 81*20f4b53cSBarry Smith monitor routine (use `NULL` if no context is desired) 82dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 83*20f4b53cSBarry Smith (may be `NULL`) 84*20f4b53cSBarry Smith 85*20f4b53cSBarry Smith Level: intermediate 86dcf2fd19SBarry Smith 87f6dfbefdSBarry Smith Note: 88dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 89f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` multiple times; all will be called in the 90dcf2fd19SBarry Smith order in which they were set. 91dcf2fd19SBarry Smith 92f6dfbefdSBarry Smith Fortran Note: 93f6dfbefdSBarry Smith Only a single monitor function can be set for each `SNESLineSearch` object 94dcf2fd19SBarry Smith 95f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()` 96dcf2fd19SBarry Smith @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **)) 98d71ae5a4SJacob Faibussowitsch { 9978064530SBarry Smith PetscInt i; 10078064530SBarry Smith PetscBool identical; 10178064530SBarry Smith 102dcf2fd19SBarry Smith PetscFunctionBegin; 103dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 10478064530SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical)); 1063ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 10778064530SBarry Smith } 1085f80ce2aSJacob Faibussowitsch PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 109dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 110dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 111dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void *)mctx; 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113dcf2fd19SBarry Smith } 114dcf2fd19SBarry Smith 115dcf2fd19SBarry Smith /*@C 116f6dfbefdSBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries 117dcf2fd19SBarry Smith 118c3339decSBarry Smith Collective 119dcf2fd19SBarry Smith 120dcf2fd19SBarry Smith Input Parameters: 121f6dfbefdSBarry Smith + ls - the `SNES` linesearch object 122f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat` 123dcf2fd19SBarry Smith 124dcf2fd19SBarry Smith Level: intermediate 125dcf2fd19SBarry Smith 126f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()` 127dcf2fd19SBarry Smith @*/ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) 129d71ae5a4SJacob Faibussowitsch { 130d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 131dcf2fd19SBarry Smith Vec Y, W, G; 132dcf2fd19SBarry Smith 133dcf2fd19SBarry Smith PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n")); 1379566063dSJacob Faibussowitsch PetscCall(VecView(Y, viewer)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n")); 1399566063dSJacob Faibussowitsch PetscCall(VecView(W, viewer)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n")); 1419566063dSJacob Faibussowitsch PetscCall(VecView(G, viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144dcf2fd19SBarry Smith } 145dcf2fd19SBarry Smith 146f40b411bSPeter Brune /*@ 147cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 148f40b411bSPeter Brune 149f6dfbefdSBarry Smith Logically Collective 150f40b411bSPeter Brune 151f40b411bSPeter Brune Input Parameters: 152f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context). 153f40b411bSPeter Brune 154f40b411bSPeter Brune Output Parameters: 1558e557f58SPeter Brune . outlinesearch - the new linesearch context 156f40b411bSPeter Brune 157162e0bf5SPeter Brune Level: developer 158162e0bf5SPeter Brune 159f6dfbefdSBarry Smith Note: 160f6dfbefdSBarry Smith The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance 161f6dfbefdSBarry Smith already associated with the `SNES`. 162f40b411bSPeter Brune 163f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()` 164f40b411bSPeter Brune @*/ 165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 166d71ae5a4SJacob Faibussowitsch { 167f1c6b773SPeter Brune SNESLineSearch linesearch; 168bf388a1fSBarry Smith 169bf7f4e0aSPeter Brune PetscFunctionBegin; 170ea5d4fccSPeter Brune PetscValidPointer(outlinesearch, 2); 1719566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 1720298fd71SBarry Smith *outlinesearch = NULL; 173f5af7f23SKarl Rupp 1749566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView)); 175bf7f4e0aSPeter Brune 1760298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1770298fd71SBarry Smith linesearch->vec_func_new = NULL; 1780298fd71SBarry Smith linesearch->vec_sol = NULL; 1790298fd71SBarry Smith linesearch->vec_func = NULL; 1800298fd71SBarry Smith linesearch->vec_update = NULL; 1819bd66eb0SPeter Brune 182bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 183bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 184bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 185bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 186422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 187bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 188bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 189bf7f4e0aSPeter Brune linesearch->damping = 1.0; 190bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 191bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 192516fe3c3SPeter Brune linesearch->rtol = 1e-8; 193516fe3c3SPeter Brune linesearch->atol = 1e-15; 194516fe3c3SPeter Brune linesearch->ltol = 1e-8; 1950298fd71SBarry Smith linesearch->precheckctx = NULL; 1960298fd71SBarry Smith linesearch->postcheckctx = NULL; 19759405d5eSPeter Brune linesearch->max_its = 1; 198bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 1993add74b1SBarry Smith linesearch->monitor = NULL; 200bf7f4e0aSPeter Brune *outlinesearch = linesearch; 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 202bf7f4e0aSPeter Brune } 203bf7f4e0aSPeter Brune 204f40b411bSPeter Brune /*@ 20578bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 20678bcb3b5SPeter Brune any required vectors. 207f40b411bSPeter Brune 208c3339decSBarry Smith Collective 209f40b411bSPeter Brune 210f40b411bSPeter Brune Input Parameters: 211f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 212f40b411bSPeter Brune 213*20f4b53cSBarry Smith Level: advanced 214*20f4b53cSBarry Smith 215f6dfbefdSBarry Smith Note: 216f6dfbefdSBarry Smith For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`. 217cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 21878bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 219f6dfbefdSBarry Smith of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be 220cd7522eaSPeter Brune allocated upfront. 221cd7522eaSPeter Brune 222f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()` 223f40b411bSPeter Brune @*/ 224f40b411bSPeter Brune 225d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 226d71ae5a4SJacob Faibussowitsch { 227bf7f4e0aSPeter Brune PetscFunctionBegin; 22848a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 229bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 23048a46eb9SPierre Jolivet if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 23148a46eb9SPierre Jolivet if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 232dbbe0bcdSBarry Smith if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup); 2339566063dSJacob Faibussowitsch if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction)); 234bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 235bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 236bf7f4e0aSPeter Brune } 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 238bf7f4e0aSPeter Brune } 239bf7f4e0aSPeter Brune 240f40b411bSPeter Brune /*@ 241f6dfbefdSBarry Smith SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search. 242f40b411bSPeter Brune 243c3339decSBarry Smith Collective 244f40b411bSPeter Brune 245f40b411bSPeter Brune Input Parameters: 246f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 247f40b411bSPeter Brune 248*20f4b53cSBarry Smith Level: developer 249*20f4b53cSBarry Smith 250f6dfbefdSBarry Smith Note: 251f6dfbefdSBarry Smith Usually only called by `SNESReset()` 252f190f2fcSBarry Smith 253f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 254f40b411bSPeter Brune @*/ 255f40b411bSPeter Brune 256d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 257d71ae5a4SJacob Faibussowitsch { 258bf7f4e0aSPeter Brune PetscFunctionBegin; 259dbbe0bcdSBarry Smith if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset); 260f5af7f23SKarl Rupp 2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_sol_new)); 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_func_new)); 263bf7f4e0aSPeter Brune 2649566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 265f5af7f23SKarl Rupp 266bf7f4e0aSPeter Brune linesearch->nwork = 0; 267bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269bf7f4e0aSPeter Brune } 270bf7f4e0aSPeter Brune 271ed07d7d7SPeter Brune /*@C 272f6dfbefdSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search 273f6dfbefdSBarry Smith ` 274ed07d7d7SPeter Brune Input Parameters: 275f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 276f6dfbefdSBarry Smith + func - function evaluation routine, this is usually the function provided with `SNESSetFunction()` 277ed07d7d7SPeter Brune 278ed07d7d7SPeter Brune Level: developer 279ed07d7d7SPeter Brune 280f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()` 281ed07d7d7SPeter Brune @*/ 282d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec)) 283d71ae5a4SJacob Faibussowitsch { 284ed07d7d7SPeter Brune PetscFunctionBegin; 285ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 286ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288ed07d7d7SPeter Brune } 289ed07d7d7SPeter Brune 29086d74e61SPeter Brune /*@C 291f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 292df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 293f190f2fcSBarry Smith determined the search direction. 29486d74e61SPeter Brune 295c3339decSBarry Smith Logically Collective 29686d74e61SPeter Brune 29786d74e61SPeter Brune Input Parameters: 298f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 299*20f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPreCheck()` 300*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 30186d74e61SPeter Brune 30286d74e61SPeter Brune Level: intermediate 30386d74e61SPeter Brune 304f6dfbefdSBarry Smith Note: 305f0b84518SBarry Smith Use `SNESLineSearchSetPostCheck()` to change the step after the line search. 306f0b84518SBarry Smith search is complete. 307f0b84518SBarry Smith 308f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 309f0b84518SBarry Smith 310f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 311f0b84518SBarry Smith `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError() 312f0b84518SBarry Smith 31386d74e61SPeter Brune @*/ 314d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx) 315d71ae5a4SJacob Faibussowitsch { 3169bd66eb0SPeter Brune PetscFunctionBegin; 317f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 318f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 31986d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32186d74e61SPeter Brune } 32286d74e61SPeter Brune 32386d74e61SPeter Brune /*@C 32453deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 32586d74e61SPeter Brune 326f899ff85SJose E. Roman Input Parameter: 327f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 32886d74e61SPeter Brune 32986d74e61SPeter Brune Output Parameters: 330*20f4b53cSBarry Smith + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchPreCheck` 331*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 33286d74e61SPeter Brune 33386d74e61SPeter Brune Level: intermediate 33486d74e61SPeter Brune 335f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 33686d74e61SPeter Brune @*/ 337d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) 338d71ae5a4SJacob Faibussowitsch { 3399bd66eb0SPeter Brune PetscFunctionBegin; 340f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 341f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 34286d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 3433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34486d74e61SPeter Brune } 34586d74e61SPeter Brune 34686d74e61SPeter Brune /*@C 347f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 348f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 34986d74e61SPeter Brune 350*20f4b53cSBarry Smith Logically Collective 35186d74e61SPeter Brune 35286d74e61SPeter Brune Input Parameters: 353f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 354*20f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPostCheck` 355*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 35686d74e61SPeter Brune 35786d74e61SPeter Brune Level: intermediate 35886d74e61SPeter Brune 359f0b84518SBarry Smith Notes: 360f0b84518SBarry Smith Use `SNESLineSearchSetPreCheck()` to change the step before the line search. 361f0b84518SBarry Smith search is complete. 362f0b84518SBarry Smith 363f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 364f0b84518SBarry Smith 365f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`, 366f0b84518SBarry Smith `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError() 36786d74e61SPeter Brune @*/ 368d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 369d71ae5a4SJacob Faibussowitsch { 37086d74e61SPeter Brune PetscFunctionBegin; 371f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 372f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 37386d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37586d74e61SPeter Brune } 37686d74e61SPeter Brune 37786d74e61SPeter Brune /*@C 378f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 37986d74e61SPeter Brune 380f899ff85SJose E. Roman Input Parameter: 381f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 38286d74e61SPeter Brune 38386d74e61SPeter Brune Output Parameters: 384f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck` 385*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 38686d74e61SPeter Brune 38786d74e61SPeter Brune Level: intermediate 38886d74e61SPeter Brune 389f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 39086d74e61SPeter Brune @*/ 391d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 392d71ae5a4SJacob Faibussowitsch { 3939bd66eb0SPeter Brune PetscFunctionBegin; 394f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 395f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 39686d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39886d74e61SPeter Brune } 39986d74e61SPeter Brune 400f40b411bSPeter Brune /*@ 401f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 402f40b411bSPeter Brune 403c3339decSBarry Smith Logically Collective 404f40b411bSPeter Brune 405f40b411bSPeter Brune Input Parameters: 4067b1df9c1SPeter Brune + linesearch - The linesearch instance. 4077b1df9c1SPeter Brune . X - The current solution 4087b1df9c1SPeter Brune - Y - The step direction 409f40b411bSPeter Brune 410f40b411bSPeter Brune Output Parameters: 4118e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 412f40b411bSPeter Brune 413f0b84518SBarry Smith Level: advanced 414f40b411bSPeter Brune 415f6dfbefdSBarry Smith Note: 416f6dfbefdSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` 417f6dfbefdSBarry Smith 418f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, 419f0b84518SBarry Smith `SNESLineSearchGetPostCheck()`` 420f40b411bSPeter Brune @*/ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) 422d71ae5a4SJacob Faibussowitsch { 423bf7f4e0aSPeter Brune PetscFunctionBegin; 424bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 4256b2b7091SBarry Smith if (linesearch->ops->precheck) { 426dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx); 42738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed, 4); 428bf7f4e0aSPeter Brune } 4293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 430bf7f4e0aSPeter Brune } 431bf7f4e0aSPeter Brune 432f40b411bSPeter Brune /*@ 433ef46b1a6SStefano Zampini SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 434f40b411bSPeter Brune 435c3339decSBarry Smith Logically Collective 436f40b411bSPeter Brune 437f40b411bSPeter Brune Input Parameters: 4387b1df9c1SPeter Brune + linesearch - The linesearch context 4397b1df9c1SPeter Brune . X - The last solution 4407b1df9c1SPeter Brune . Y - The step direction 4417b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 442f40b411bSPeter Brune 443f40b411bSPeter Brune Output Parameters: 44478bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 44578bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 446f40b411bSPeter Brune 447*20f4b53cSBarry Smith Level: developer 448*20f4b53cSBarry Smith 449f6dfbefdSBarry Smith Note: 450f6dfbefdSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` 451f6dfbefdSBarry Smith 452db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 453f40b411bSPeter Brune @*/ 454d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 455d71ae5a4SJacob Faibussowitsch { 456bf7f4e0aSPeter Brune PetscFunctionBegin; 457bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 458bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4596b2b7091SBarry Smith if (linesearch->ops->postcheck) { 460dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx); 46138bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5); 46238bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6); 46386d74e61SPeter Brune } 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46586d74e61SPeter Brune } 46686d74e61SPeter Brune 46786d74e61SPeter Brune /*@C 46886d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 46986d74e61SPeter Brune 470c3339decSBarry Smith Logically Collective 47186d74e61SPeter Brune 4724165533cSJose E. Roman Input Parameters: 47386d74e61SPeter Brune + linesearch - linesearch context 47486d74e61SPeter Brune . X - base state for this step 475907376e6SBarry Smith - ctx - context for this function 47686d74e61SPeter Brune 47797bb3fdcSJose E. Roman Input/Output Parameter: 47897bb3fdcSJose E. Roman . Y - correction, possibly modified 47997bb3fdcSJose E. Roman 48097bb3fdcSJose E. Roman Output Parameter: 48197bb3fdcSJose E. Roman . changed - flag indicating that Y was modified 48286d74e61SPeter Brune 48386d74e61SPeter Brune Options Database Key: 484cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 485cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 48686d74e61SPeter Brune 48786d74e61SPeter Brune Level: advanced 48886d74e61SPeter Brune 48986d74e61SPeter Brune Notes: 490f6dfbefdSBarry Smith This function should be passed to `SNESLineSearchSetPreCheck()` 49186d74e61SPeter Brune 49286d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 49386d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 49486d74e61SPeter Brune is generally not useful when using a Newton linearization. 49586d74e61SPeter Brune 49686d74e61SPeter Brune Reference: 497f6dfbefdSBarry Smith . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 49886d74e61SPeter Brune 499f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()` 50086d74e61SPeter Brune @*/ 501d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) 502d71ae5a4SJacob Faibussowitsch { 50386d74e61SPeter Brune PetscReal angle = *(PetscReal *)linesearch->precheckctx; 50486d74e61SPeter Brune Vec Ylast; 50586d74e61SPeter Brune PetscScalar dot; 50686d74e61SPeter Brune PetscInt iter; 50786d74e61SPeter Brune PetscReal ynorm, ylastnorm, theta, angle_radians; 50886d74e61SPeter Brune SNES snes; 50986d74e61SPeter Brune 51086d74e61SPeter Brune PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 5129566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast)); 51386d74e61SPeter Brune if (!Ylast) { 5149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Ylast)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast)); 5169566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Ylast)); 51786d74e61SPeter Brune } 5189566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 51986d74e61SPeter Brune if (iter < 2) { 5209566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 52186d74e61SPeter Brune *changed = PETSC_FALSE; 5223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52386d74e61SPeter Brune } 52486d74e61SPeter Brune 5259566063dSJacob Faibussowitsch PetscCall(VecDot(Y, Ylast, &dot)); 5269566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); 5279566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 52886d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 529255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0)); 53086d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 53186d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 53286d74e61SPeter Brune /* Modify the step Y */ 53386d74e61SPeter Brune PetscReal alpha, ydiffnorm; 5349566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ylast, -1.0, Y)); 5359566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm)); 536f85e2ce2SBarry Smith alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0; 5379566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 5389566063dSJacob Faibussowitsch PetscCall(VecScale(Y, alpha)); 5399566063dSJacob Faibussowitsch PetscCall(PetscInfo(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)); 540a47ec85fSBarry Smith *changed = PETSC_TRUE; 54186d74e61SPeter Brune } else { 5429566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle)); 5439566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 54486d74e61SPeter Brune *changed = PETSC_FALSE; 545bf7f4e0aSPeter Brune } 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 547bf7f4e0aSPeter Brune } 548bf7f4e0aSPeter Brune 549f40b411bSPeter Brune /*@ 550cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 551f40b411bSPeter Brune 552c3339decSBarry Smith Collective 553f40b411bSPeter Brune 554f40b411bSPeter Brune Input Parameters: 5558e557f58SPeter Brune + linesearch - The linesearch context 5568e557f58SPeter Brune - Y - The search direction 557f40b411bSPeter Brune 5586b867d5aSJose E. Roman Input/Output Parameters: 5596b867d5aSJose E. Roman + X - The current solution, on output the new solution 5606b867d5aSJose E. Roman . F - The current function, on output the new function 5616b867d5aSJose E. Roman - fnorm - The current norm, on output the new function norm 562f40b411bSPeter Brune 563cd7522eaSPeter Brune Options Database Keys: 5640b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell 565dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 5661fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 5671fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 5683c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 5693c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 570cd7522eaSPeter Brune 571*20f4b53cSBarry Smith Level: Intermediate 572*20f4b53cSBarry Smith 573cd7522eaSPeter Brune Notes: 574f6dfbefdSBarry Smith This is typically called from within a `SNESSolve()` implementation in order to 575f6dfbefdSBarry Smith help with convergence of the nonlinear method. Various `SNES` types use line searches 576cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 577cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 578f6dfbefdSBarry Smith application of the line search may invoke `SNESComputeFunction()` several times, and 579cd7522eaSPeter Brune therefore may be fairly expensive. 580cd7522eaSPeter Brune 581f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`, 582db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetType()` 583f40b411bSPeter Brune @*/ 584d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) 585d71ae5a4SJacob Faibussowitsch { 586bf388a1fSBarry Smith PetscFunctionBegin; 587f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 588bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 589bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 590064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Y, VEC_CLASSID, 5); 591bf7f4e0aSPeter Brune 592422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 593bf7f4e0aSPeter Brune 594bf7f4e0aSPeter Brune linesearch->vec_sol = X; 595bf7f4e0aSPeter Brune linesearch->vec_update = Y; 596bf7f4e0aSPeter Brune linesearch->vec_func = F; 597bf7f4e0aSPeter Brune 5989566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(linesearch)); 599bf7f4e0aSPeter Brune 600f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 601bf7f4e0aSPeter Brune 602f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 60348a46eb9SPierre Jolivet else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 604bf7f4e0aSPeter Brune 6059566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 606bf7f4e0aSPeter Brune 607dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, apply); 608bf7f4e0aSPeter Brune 6099566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 610bf7f4e0aSPeter Brune 611f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 613bf7f4e0aSPeter Brune } 614bf7f4e0aSPeter Brune 615f40b411bSPeter Brune /*@ 616f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 617f40b411bSPeter Brune 618c3339decSBarry Smith Collective 619f40b411bSPeter Brune 620f40b411bSPeter Brune Input Parameters: 6218e557f58SPeter Brune . linesearch - The linesearch context 622f40b411bSPeter Brune 62384238204SBarry Smith Level: developer 624f40b411bSPeter Brune 625f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()` 626f40b411bSPeter Brune @*/ 627d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) 628d71ae5a4SJacob Faibussowitsch { 629bf7f4e0aSPeter Brune PetscFunctionBegin; 6303ba16761SJacob Faibussowitsch if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS); 631f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1); 6329371c9d4SSatish Balay if (--((PetscObject)(*linesearch))->refct > 0) { 6339371c9d4SSatish Balay *linesearch = NULL; 6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6359371c9d4SSatish Balay } 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch)); 6379566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset(*linesearch)); 6383ba16761SJacob Faibussowitsch PetscTryTypeMethod(*linesearch, destroy); 6399566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 6409566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorCancel((*linesearch))); 6419566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(linesearch)); 6423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 643bf7f4e0aSPeter Brune } 644bf7f4e0aSPeter Brune 645f40b411bSPeter Brune /*@ 646dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 647bf7f4e0aSPeter Brune 648c3339decSBarry Smith Logically Collective 649f6dfbefdSBarry Smith 650bf7f4e0aSPeter Brune Input Parameters: 651dcf2fd19SBarry Smith + linesearch - the linesearch object 652*20f4b53cSBarry Smith - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor 653bf7f4e0aSPeter Brune 654f6dfbefdSBarry Smith Options Database Key: 655dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 656bf7f4e0aSPeter Brune 657bf7f4e0aSPeter Brune Level: intermediate 658bf7f4e0aSPeter Brune 659f6dfbefdSBarry Smith Developer Note: 660f6dfbefdSBarry Smith This monitor is implemented differently than the other line search monitors that are set with 661f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the 662d12e167eSBarry Smith line search that are not visible to the other monitors. 663dcf2fd19SBarry Smith 664f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`, 665f6dfbefdSBarry Smith `SNESLineSearchMonitorSetFromOptions()` 666bf7f4e0aSPeter Brune @*/ 667d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 668d71ae5a4SJacob Faibussowitsch { 669bf7f4e0aSPeter Brune PetscFunctionBegin; 6709566063dSJacob Faibussowitsch if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer)); 6719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&linesearch->monitor)); 672dcf2fd19SBarry Smith linesearch->monitor = viewer; 6733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 674bf7f4e0aSPeter Brune } 675bf7f4e0aSPeter Brune 676f40b411bSPeter Brune /*@ 677f6dfbefdSBarry Smith SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor. 678f6dfbefdSBarry Smith 679c3339decSBarry Smith Logically Collective 6806a388c36SPeter Brune 681f190f2fcSBarry Smith Input Parameter: 6828e557f58SPeter Brune . linesearch - linesearch context 683f40b411bSPeter Brune 684f190f2fcSBarry Smith Output Parameter: 6858e557f58SPeter Brune . monitor - monitor context 686f40b411bSPeter Brune 687f40b411bSPeter Brune Level: intermediate 688f40b411bSPeter Brune 689f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 690f40b411bSPeter Brune @*/ 691d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 692d71ae5a4SJacob Faibussowitsch { 6936a388c36SPeter Brune PetscFunctionBegin; 694f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 6956a388c36SPeter Brune *monitor = linesearch->monitor; 6963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6976a388c36SPeter Brune } 6986a388c36SPeter Brune 699dcf2fd19SBarry Smith /*@C 700dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 701dcf2fd19SBarry Smith 702c3339decSBarry Smith Collective 703dcf2fd19SBarry Smith 704dcf2fd19SBarry Smith Input Parameters: 705f6dfbefdSBarry Smith + ls - `SNESLineSearch` object you wish to monitor 706dcf2fd19SBarry Smith . name - the monitor type one is seeking 707dcf2fd19SBarry Smith . help - message indicating what monitoring is done 708dcf2fd19SBarry Smith . manual - manual page for the monitor 709dcf2fd19SBarry Smith . monitor - the monitor function 710f6dfbefdSBarry 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` 711dcf2fd19SBarry Smith 712dcf2fd19SBarry Smith Level: developer 713dcf2fd19SBarry Smith 714f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 715db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 716db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 717db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 718c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 719db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 720db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 721dcf2fd19SBarry Smith @*/ 722d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *)) 723d71ae5a4SJacob Faibussowitsch { 724dcf2fd19SBarry Smith PetscViewer viewer; 725dcf2fd19SBarry Smith PetscViewerFormat format; 726dcf2fd19SBarry Smith PetscBool flg; 727dcf2fd19SBarry Smith 728dcf2fd19SBarry Smith PetscFunctionBegin; 7299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg)); 730dcf2fd19SBarry Smith if (flg) { 731d12e167eSBarry Smith PetscViewerAndFormat *vf; 7329566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 7339566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 7341baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls, vf)); 7359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 736dcf2fd19SBarry Smith } 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 738dcf2fd19SBarry Smith } 739dcf2fd19SBarry Smith 740f40b411bSPeter Brune /*@ 741f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 742f40b411bSPeter Brune 743c3339decSBarry Smith Logically Collective 744f6dfbefdSBarry Smith 745f899ff85SJose E. Roman Input Parameter: 7468e557f58SPeter Brune . linesearch - linesearch context 747f40b411bSPeter Brune 748f40b411bSPeter Brune Options Database Keys: 7490b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 7503c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 751f6dfbefdSBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`) 75271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 7531a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 7541a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 7551a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 7561a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 7571a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 758dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 759dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 7608e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 761cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 7621a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 763d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 764f40b411bSPeter Brune 765f40b411bSPeter Brune Level: intermediate 766f40b411bSPeter Brune 767f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 768db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 769f40b411bSPeter Brune @*/ 770d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 771d71ae5a4SJacob Faibussowitsch { 7721a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 773bf7f4e0aSPeter Brune char type[256]; 774bf7f4e0aSPeter Brune PetscBool flg, set; 775dcf2fd19SBarry Smith PetscViewer viewer; 776bf388a1fSBarry Smith 777bf7f4e0aSPeter Brune PetscFunctionBegin; 7789566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll()); 779bf7f4e0aSPeter Brune 780d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch); 781f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 7829566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg)); 783bf7f4e0aSPeter Brune if (flg) { 7849566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, type)); 785bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 7869566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, deft)); 787bf7f4e0aSPeter Brune } 788bf7f4e0aSPeter Brune 7899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set)); 790dcf2fd19SBarry Smith if (set) { 7919566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer)); 7929566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 793dcf2fd19SBarry Smith } 7949566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL)); 795bf7f4e0aSPeter Brune 7961a9b3a06SPeter Brune /* tolerances */ 7979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL)); 7989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL)); 7999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL)); 8009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL)); 8019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL)); 8029566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL)); 8037a35526eSPeter Brune 8041a9b3a06SPeter Brune /* damping parameters */ 8059566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 8061a9b3a06SPeter Brune 8079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 8081a9b3a06SPeter Brune 8091a9b3a06SPeter Brune /* precheck */ 8109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 8117a35526eSPeter Brune if (set) { 8127a35526eSPeter Brune if (flg) { 8137a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 814f5af7f23SKarl Rupp 815d0609cedSBarry Smith PetscCall(PetscOptionsReal("-snes_linesearch_precheck_picard_angle", "Maximum angle at which to activate the correction", "none", linesearch->precheck_picard_angle, &linesearch->precheck_picard_angle, NULL)); 8169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 8177a35526eSPeter Brune } else { 8189566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 8197a35526eSPeter Brune } 8207a35526eSPeter Brune } 8219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 8229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 8237a35526eSPeter Brune 824dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 8255a9b6599SPeter Brune 826dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 827d0609cedSBarry Smith PetscOptionsEnd(); 8283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 829bf7f4e0aSPeter Brune } 830bf7f4e0aSPeter Brune 831f40b411bSPeter Brune /*@ 832f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 833f40b411bSPeter Brune 834*20f4b53cSBarry Smith Logically Collective 835*20f4b53cSBarry Smith 836f40b411bSPeter Brune Input Parameters: 8378e557f58SPeter Brune . linesearch - linesearch context 838f40b411bSPeter Brune 839f40b411bSPeter Brune Level: intermediate 840f40b411bSPeter Brune 841f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 842f40b411bSPeter Brune @*/ 843d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 844d71ae5a4SJacob Faibussowitsch { 8457f1410a3SPeter Brune PetscBool iascii; 846bf388a1fSBarry Smith 847bf7f4e0aSPeter Brune PetscFunctionBegin; 8487f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 84948a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 8507f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 8517f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2); 852f40b411bSPeter Brune 8539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 8547f1410a3SPeter Brune if (iascii) { 8559566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 8569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 857dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer); 8589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 8599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 8609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 86163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 8626b2b7091SBarry Smith if (linesearch->ops->precheck) { 8636b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 86463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 8657f1410a3SPeter Brune } else { 86663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 8677f1410a3SPeter Brune } 8687f1410a3SPeter Brune } 86948a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 8707f1410a3SPeter Brune } 8713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 872bf7f4e0aSPeter Brune } 873bf7f4e0aSPeter Brune 874ea5d4fccSPeter Brune /*@C 875a80ff896SJed Brown SNESLineSearchGetType - Gets the linesearch type 876a80ff896SJed Brown 877c3339decSBarry Smith Logically Collective 878a80ff896SJed Brown 879a80ff896SJed Brown Input Parameters: 880a80ff896SJed Brown . linesearch - linesearch context 881a80ff896SJed Brown 882a80ff896SJed Brown Output Parameters: 883*20f4b53cSBarry Smith - type - The type of line search, or `NULL` if not set 884a80ff896SJed Brown 885a80ff896SJed Brown Level: intermediate 886a80ff896SJed Brown 887f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 888a80ff896SJed Brown @*/ 889d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 890d71ae5a4SJacob Faibussowitsch { 891a80ff896SJed Brown PetscFunctionBegin; 892a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 893dadcf809SJacob Faibussowitsch PetscValidPointer(type, 2); 894a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 8953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 896a80ff896SJed Brown } 897a80ff896SJed Brown 898a80ff896SJed Brown /*@C 899f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 900f40b411bSPeter Brune 901c3339decSBarry Smith Logically Collective 902f190f2fcSBarry Smith 903f40b411bSPeter Brune Input Parameters: 9048e557f58SPeter Brune + linesearch - linesearch context 905f40b411bSPeter Brune - type - The type of line search to be used 906f40b411bSPeter Brune 907cd7522eaSPeter Brune Available Types: 908f6dfbefdSBarry Smith + `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step 909f6dfbefdSBarry Smith . `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function 910f6dfbefdSBarry Smith . `SNESLINESEARCHL2` - Secant line search over the L2 norm of the function 911f6dfbefdSBarry Smith . `SNESLINESEARCHCP` - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 912f6dfbefdSBarry Smith . `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch 913f6dfbefdSBarry Smith - `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation 9141fe24845SBarry Smith 9153c7db156SBarry Smith Options Database Key: 9160b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 917cd7522eaSPeter Brune 918f40b411bSPeter Brune Level: intermediate 919f40b411bSPeter Brune 920f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 921f40b411bSPeter Brune @*/ 922d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 923d71ae5a4SJacob Faibussowitsch { 924bf7f4e0aSPeter Brune PetscBool match; 9255f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 926bf7f4e0aSPeter Brune 927bf7f4e0aSPeter Brune PetscFunctionBegin; 928f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 929bf7f4e0aSPeter Brune PetscValidCharPointer(type, 2); 930bf7f4e0aSPeter Brune 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 9323ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 933bf7f4e0aSPeter Brune 9349566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 9355f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 936bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 937bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 9389566063dSJacob Faibussowitsch PetscCall((*(linesearch)->ops->destroy)(linesearch)); 9390298fd71SBarry Smith linesearch->ops->destroy = NULL; 940bf7f4e0aSPeter Brune } 941f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9429e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9439e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9449e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9459e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 946bf7f4e0aSPeter Brune 9479566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 9489566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 9493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 950bf7f4e0aSPeter Brune } 951bf7f4e0aSPeter Brune 952f40b411bSPeter Brune /*@ 953f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 954f40b411bSPeter Brune 955f40b411bSPeter Brune Input Parameters: 9568e557f58SPeter Brune + linesearch - linesearch context 957f40b411bSPeter Brune - snes - The snes instance 958f40b411bSPeter Brune 95978bcb3b5SPeter Brune Level: developer 96078bcb3b5SPeter Brune 961f6dfbefdSBarry Smith Note: 962f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 963f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 96478bcb3b5SPeter Brune implementations. 965f40b411bSPeter Brune 966f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 967f40b411bSPeter Brune @*/ 968d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 969d71ae5a4SJacob Faibussowitsch { 970bf7f4e0aSPeter Brune PetscFunctionBegin; 971f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 972bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 973bf7f4e0aSPeter Brune linesearch->snes = snes; 9743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 975bf7f4e0aSPeter Brune } 976bf7f4e0aSPeter Brune 977f40b411bSPeter Brune /*@ 978f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 979f6dfbefdSBarry Smith Having an associated `SNES` is necessary because most line search implementations must be able to 980f6dfbefdSBarry Smith evaluate the function using `SNESComputeFunction()` for the associated `SNES`. This routine 981f6dfbefdSBarry Smith is used in the line search implementations when one must get this associated `SNES` instance. 982f6dfbefdSBarry Smith 983f6dfbefdSBarry Smith Not Collective 984f40b411bSPeter Brune 985f40b411bSPeter Brune Input Parameters: 9868e557f58SPeter Brune . linesearch - linesearch context 987f40b411bSPeter Brune 988f40b411bSPeter Brune Output Parameters: 989f40b411bSPeter Brune . snes - The snes instance 990f40b411bSPeter Brune 9918141a3b9SPeter Brune Level: developer 992f40b411bSPeter Brune 993f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 994f40b411bSPeter Brune @*/ 995d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 996d71ae5a4SJacob Faibussowitsch { 997bf7f4e0aSPeter Brune PetscFunctionBegin; 998f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9996a388c36SPeter Brune PetscValidPointer(snes, 2); 1000bf7f4e0aSPeter Brune *snes = linesearch->snes; 10013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1002bf7f4e0aSPeter Brune } 1003bf7f4e0aSPeter Brune 1004f40b411bSPeter Brune /*@ 1005f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1006f40b411bSPeter Brune 1007f6dfbefdSBarry Smith Not Collective 1008f6dfbefdSBarry Smith 1009f40b411bSPeter Brune Input Parameters: 10108e557f58SPeter Brune . linesearch - linesearch context 1011f40b411bSPeter Brune 1012f40b411bSPeter Brune Output Parameters: 1013f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()` 1014f40b411bSPeter Brune 101578bcb3b5SPeter Brune Level: advanced 101678bcb3b5SPeter Brune 1017f6dfbefdSBarry Smith Note: 10188e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 101978bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 1020f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the 102178bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 102278bcb3b5SPeter Brune 1023f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1024f40b411bSPeter Brune @*/ 1025d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1026d71ae5a4SJacob Faibussowitsch { 10276a388c36SPeter Brune PetscFunctionBegin; 1028f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1029534a8f05SLisandro Dalcin PetscValidRealPointer(lambda, 2); 10306a388c36SPeter Brune *lambda = linesearch->lambda; 10313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10326a388c36SPeter Brune } 10336a388c36SPeter Brune 1034f40b411bSPeter Brune /*@ 1035f6dfbefdSBarry Smith SNESLineSearchSetLambda - Sets the linesearch steplength 1036f40b411bSPeter Brune 1037f40b411bSPeter Brune Input Parameters: 10388e557f58SPeter Brune + linesearch - linesearch context 1039f40b411bSPeter Brune - lambda - The last steplength. 1040f40b411bSPeter Brune 1041*20f4b53cSBarry Smith Level: advanced 1042*20f4b53cSBarry Smith 1043f6dfbefdSBarry Smith Note: 1044f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()` 1045f6dfbefdSBarry Smith to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1046cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1047cd7522eaSPeter Brune as an inner scaling parameter. 1048cd7522eaSPeter Brune 1049f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()` 1050f40b411bSPeter Brune @*/ 1051d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1052d71ae5a4SJacob Faibussowitsch { 10536a388c36SPeter Brune PetscFunctionBegin; 1054f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10556a388c36SPeter Brune linesearch->lambda = lambda; 10563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10576a388c36SPeter Brune } 10586a388c36SPeter Brune 1059f40b411bSPeter Brune /*@ 10603c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 106178bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 106278bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 106378bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1064f40b411bSPeter Brune 1065f6dfbefdSBarry Smith Not Collective 1066f6dfbefdSBarry Smith 1067f899ff85SJose E. Roman Input Parameter: 10688e557f58SPeter Brune . linesearch - linesearch context 1069f40b411bSPeter Brune 1070f40b411bSPeter Brune Output Parameters: 1071516fe3c3SPeter Brune + steptol - The minimum steplength 10726cc8e53bSPeter Brune . maxstep - The maximum steplength 1073516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1074516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1075516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1076516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1077f40b411bSPeter Brune 107878bcb3b5SPeter Brune Level: intermediate 107978bcb3b5SPeter Brune 1080f6dfbefdSBarry Smith Note: 108178bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 10823c7d6663SPeter Brune the type requires. 1083516fe3c3SPeter Brune 1084f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1085f40b411bSPeter Brune @*/ 1086d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1087d71ae5a4SJacob Faibussowitsch { 10886a388c36SPeter Brune PetscFunctionBegin; 1089f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1090516fe3c3SPeter Brune if (steptol) { 1091534a8f05SLisandro Dalcin PetscValidRealPointer(steptol, 2); 10926a388c36SPeter Brune *steptol = linesearch->steptol; 1093516fe3c3SPeter Brune } 1094516fe3c3SPeter Brune if (maxstep) { 1095534a8f05SLisandro Dalcin PetscValidRealPointer(maxstep, 3); 1096516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1097516fe3c3SPeter Brune } 1098516fe3c3SPeter Brune if (rtol) { 1099534a8f05SLisandro Dalcin PetscValidRealPointer(rtol, 4); 1100516fe3c3SPeter Brune *rtol = linesearch->rtol; 1101516fe3c3SPeter Brune } 1102516fe3c3SPeter Brune if (atol) { 1103534a8f05SLisandro Dalcin PetscValidRealPointer(atol, 5); 1104516fe3c3SPeter Brune *atol = linesearch->atol; 1105516fe3c3SPeter Brune } 1106516fe3c3SPeter Brune if (ltol) { 1107534a8f05SLisandro Dalcin PetscValidRealPointer(ltol, 6); 1108516fe3c3SPeter Brune *ltol = linesearch->ltol; 1109516fe3c3SPeter Brune } 1110516fe3c3SPeter Brune if (max_its) { 1111534a8f05SLisandro Dalcin PetscValidIntPointer(max_its, 7); 1112516fe3c3SPeter Brune *max_its = linesearch->max_its; 1113516fe3c3SPeter Brune } 11143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11156a388c36SPeter Brune } 11166a388c36SPeter Brune 1117f40b411bSPeter Brune /*@ 11183c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 111978bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 112078bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 112178bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1122f40b411bSPeter Brune 1123f6dfbefdSBarry Smith Collective 1124f6dfbefdSBarry Smith 1125f40b411bSPeter Brune Input Parameters: 11268e557f58SPeter Brune + linesearch - linesearch context 1127516fe3c3SPeter Brune . steptol - The minimum steplength 11286cc8e53bSPeter Brune . maxstep - The maximum steplength 1129516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1130516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1131516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1132516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1133f40b411bSPeter Brune 1134*20f4b53cSBarry Smith Level: intermediate 1135*20f4b53cSBarry Smith 1136f6dfbefdSBarry Smith Note: 1137f6dfbefdSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in 113878bcb3b5SPeter Brune place of an argument. 1139f40b411bSPeter Brune 1140f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1141f40b411bSPeter Brune @*/ 1142d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 1143d71ae5a4SJacob Faibussowitsch { 11446a388c36SPeter Brune PetscFunctionBegin; 1145f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1146d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1147d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1148d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1149d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1150d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1151d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch, max_its, 7); 1152d3952184SSatish Balay 115313bcc0bdSJacob Faibussowitsch if (steptol != (PetscReal)PETSC_DEFAULT) { 11545f80ce2aSJacob Faibussowitsch PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 11556a388c36SPeter Brune linesearch->steptol = steptol; 1156d3952184SSatish Balay } 1157d3952184SSatish Balay 115813bcc0bdSJacob Faibussowitsch if (maxstep != (PetscReal)PETSC_DEFAULT) { 11595f80ce2aSJacob Faibussowitsch PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1160516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1161d3952184SSatish Balay } 1162d3952184SSatish Balay 116313bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) { 11642061ca32SJunchao Zhang PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %14.12e must be non-negative and less than 1.0", (double)rtol); 1165516fe3c3SPeter Brune linesearch->rtol = rtol; 1166d3952184SSatish Balay } 1167d3952184SSatish Balay 116813bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) { 11695f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1170516fe3c3SPeter Brune linesearch->atol = atol; 1171d3952184SSatish Balay } 1172d3952184SSatish Balay 117313bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) { 11745f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1175516fe3c3SPeter Brune linesearch->ltol = ltol; 1176d3952184SSatish Balay } 1177d3952184SSatish Balay 1178d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 117963a3b9bcSJacob Faibussowitsch PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its); 1180516fe3c3SPeter Brune linesearch->max_its = max_its; 1181d3952184SSatish Balay } 11823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11836a388c36SPeter Brune } 11846a388c36SPeter Brune 1185f40b411bSPeter Brune /*@ 1186f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1187f40b411bSPeter Brune 1188f40b411bSPeter Brune Input Parameters: 11898e557f58SPeter Brune . linesearch - linesearch context 1190f40b411bSPeter Brune 1191f40b411bSPeter Brune Output Parameters: 11928e557f58SPeter Brune . damping - The damping parameter 1193f40b411bSPeter Brune 119478bcb3b5SPeter Brune Level: advanced 1195f40b411bSPeter Brune 1196db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN` 1197f40b411bSPeter Brune @*/ 1198f40b411bSPeter Brune 1199d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1200d71ae5a4SJacob Faibussowitsch { 12016a388c36SPeter Brune PetscFunctionBegin; 1202f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1203534a8f05SLisandro Dalcin PetscValidRealPointer(damping, 2); 12046a388c36SPeter Brune *damping = linesearch->damping; 12053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12066a388c36SPeter Brune } 12076a388c36SPeter Brune 1208f40b411bSPeter Brune /*@ 1209fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1210f40b411bSPeter Brune 1211f40b411bSPeter Brune Input Parameters: 121203fd4120SBarry Smith + linesearch - linesearch context 121303fd4120SBarry Smith - damping - The damping parameter 1214f40b411bSPeter Brune 1215f6dfbefdSBarry Smith Options Database Key: 1216f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value 1217f6dfbefdSBarry Smith 1218f40b411bSPeter Brune Level: intermediate 1219f40b411bSPeter Brune 1220f6dfbefdSBarry Smith Note: 1221f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1222cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 122378bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1224cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1225cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1226cd7522eaSPeter Brune 1227f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()` 1228f40b411bSPeter Brune @*/ 1229d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1230d71ae5a4SJacob Faibussowitsch { 12316a388c36SPeter Brune PetscFunctionBegin; 1232f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12336a388c36SPeter Brune linesearch->damping = damping; 12343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12356a388c36SPeter Brune } 12366a388c36SPeter Brune 123759405d5eSPeter Brune /*@ 123859405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 123959405d5eSPeter Brune 1240f6dfbefdSBarry Smith Input Parameter: 124178bcb3b5SPeter Brune . linesearch - linesearch context 124259405d5eSPeter Brune 1243f6dfbefdSBarry Smith Output Parameter: 12448e557f58SPeter Brune . order - The order 124559405d5eSPeter Brune 124678bcb3b5SPeter Brune Possible Values for order: 1247f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1248f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1249f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 125078bcb3b5SPeter Brune 125159405d5eSPeter Brune Level: intermediate 125259405d5eSPeter Brune 1253f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()` 125459405d5eSPeter Brune @*/ 125559405d5eSPeter Brune 1256d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1257d71ae5a4SJacob Faibussowitsch { 125859405d5eSPeter Brune PetscFunctionBegin; 125959405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1260534a8f05SLisandro Dalcin PetscValidIntPointer(order, 2); 126159405d5eSPeter Brune *order = linesearch->order; 12623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126359405d5eSPeter Brune } 126459405d5eSPeter Brune 126559405d5eSPeter Brune /*@ 12661f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 126759405d5eSPeter Brune 126859405d5eSPeter Brune Input Parameters: 1269a2b725a8SWilliam Gropp + linesearch - linesearch context 1270a2b725a8SWilliam Gropp - order - The damping parameter 127159405d5eSPeter Brune 127259405d5eSPeter Brune Level: intermediate 127359405d5eSPeter Brune 127478bcb3b5SPeter Brune Possible Values for order: 1275f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1276f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1277f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 127878bcb3b5SPeter Brune 127959405d5eSPeter Brune Notes: 128059405d5eSPeter Brune Variable orders are supported by the following line searches: 128178bcb3b5SPeter Brune + bt - cubic and quadratic 128278bcb3b5SPeter Brune - cp - linear and quadratic 128359405d5eSPeter Brune 1284f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 128559405d5eSPeter Brune @*/ 1286d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1287d71ae5a4SJacob Faibussowitsch { 128859405d5eSPeter Brune PetscFunctionBegin; 128959405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 129059405d5eSPeter Brune linesearch->order = order; 12913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129259405d5eSPeter Brune } 129359405d5eSPeter Brune 1294f40b411bSPeter Brune /*@ 1295f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1296f40b411bSPeter Brune 1297f6dfbefdSBarry Smith Not Collective 1298f6dfbefdSBarry Smith 1299f899ff85SJose E. Roman Input Parameter: 130078bcb3b5SPeter Brune . linesearch - linesearch context 1301f40b411bSPeter Brune 1302f40b411bSPeter Brune Output Parameters: 1303f40b411bSPeter Brune + xnorm - The norm of the current solution 1304f40b411bSPeter Brune . fnorm - The norm of the current function 1305f40b411bSPeter Brune - ynorm - The norm of the current update 1306f40b411bSPeter Brune 130778bcb3b5SPeter Brune Level: developer 1308f40b411bSPeter Brune 1309f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1310f40b411bSPeter Brune @*/ 1311d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1312d71ae5a4SJacob Faibussowitsch { 1313bf7f4e0aSPeter Brune PetscFunctionBegin; 1314f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1315f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1316f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1317f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 13183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1319bf7f4e0aSPeter Brune } 1320bf7f4e0aSPeter Brune 1321f40b411bSPeter Brune /*@ 1322f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1323f40b411bSPeter Brune 1324c3339decSBarry Smith Collective 1325f6dfbefdSBarry Smith 1326f40b411bSPeter Brune Input Parameters: 132778bcb3b5SPeter Brune + linesearch - linesearch context 1328f40b411bSPeter Brune . xnorm - The norm of the current solution 1329f40b411bSPeter Brune . fnorm - The norm of the current function 1330f40b411bSPeter Brune - ynorm - The norm of the current update 1331f40b411bSPeter Brune 1332f6dfbefdSBarry Smith Level: developer 1333f40b411bSPeter Brune 1334f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1335f40b411bSPeter Brune @*/ 1336d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1337d71ae5a4SJacob Faibussowitsch { 13386a388c36SPeter Brune PetscFunctionBegin; 1339f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 13406a388c36SPeter Brune linesearch->xnorm = xnorm; 13416a388c36SPeter Brune linesearch->fnorm = fnorm; 13426a388c36SPeter Brune linesearch->ynorm = ynorm; 13433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13446a388c36SPeter Brune } 13456a388c36SPeter Brune 1346f40b411bSPeter Brune /*@ 1347f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1348f40b411bSPeter Brune 1349f40b411bSPeter Brune Input Parameters: 135078bcb3b5SPeter Brune . linesearch - linesearch context 1351f40b411bSPeter Brune 1352*20f4b53cSBarry Smith Options Database Key: 13538e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1354f40b411bSPeter Brune 1355f40b411bSPeter Brune Level: intermediate 1356f40b411bSPeter Brune 1357f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1358f40b411bSPeter Brune @*/ 1359d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1360d71ae5a4SJacob Faibussowitsch { 13619bd66eb0SPeter Brune SNES snes; 1362bf388a1fSBarry Smith 13636a388c36SPeter Brune PetscFunctionBegin; 13646a388c36SPeter Brune if (linesearch->norms) { 13659bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 13669566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 13679566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13689566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13699566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 13709bd66eb0SPeter Brune } else { 13719566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13729566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13739566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13749566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13759566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13769566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13776a388c36SPeter Brune } 13789bd66eb0SPeter Brune } 13793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13806a388c36SPeter Brune } 13816a388c36SPeter Brune 13826f263ca3SPeter Brune /*@ 13836f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 13846f263ca3SPeter Brune 13856f263ca3SPeter Brune Input Parameters: 138678bcb3b5SPeter Brune + linesearch - linesearch context 138778bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 13886f263ca3SPeter Brune 1389*20f4b53cSBarry Smith Options Database Key: 13900b00b554SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch 13916f263ca3SPeter Brune 1392*20f4b53cSBarry Smith Level: intermediate 1393*20f4b53cSBarry Smith 1394f6dfbefdSBarry Smith Note: 1395f6dfbefdSBarry Smith This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm. 13966f263ca3SPeter Brune 1397f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 13986f263ca3SPeter Brune @*/ 1399d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1400d71ae5a4SJacob Faibussowitsch { 14016f263ca3SPeter Brune PetscFunctionBegin; 14026f263ca3SPeter Brune linesearch->norms = flg; 14033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14046f263ca3SPeter Brune } 14056f263ca3SPeter Brune 1406f40b411bSPeter Brune /*@ 1407f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1408f6dfbefdSBarry Smith 1409f6dfbefdSBarry Smith Not Collective on but the vectors are parallel 1410f40b411bSPeter Brune 1411f899ff85SJose E. Roman Input Parameter: 141278bcb3b5SPeter Brune . linesearch - linesearch context 1413f40b411bSPeter Brune 1414f40b411bSPeter Brune Output Parameters: 14156232e825SPeter Brune + X - Solution vector 14166232e825SPeter Brune . F - Function vector 14176232e825SPeter Brune . Y - Search direction vector 14186232e825SPeter Brune . W - Solution work vector 14196232e825SPeter Brune - G - Function work vector 14206232e825SPeter Brune 1421*20f4b53cSBarry Smith Level: advanced 1422*20f4b53cSBarry Smith 14237bba9028SPeter Brune Notes: 1424*20f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a 1425*20f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the 1426*20f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the 14276232e825SPeter Brune function evaluated at the new solution. 1428f40b411bSPeter Brune 1429f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 14302a7a6963SBarry Smith 1431f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1432f40b411bSPeter Brune @*/ 1433d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1434d71ae5a4SJacob Faibussowitsch { 14356a388c36SPeter Brune PetscFunctionBegin; 1436f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14376a388c36SPeter Brune if (X) { 14386a388c36SPeter Brune PetscValidPointer(X, 2); 14396a388c36SPeter Brune *X = linesearch->vec_sol; 14406a388c36SPeter Brune } 14416a388c36SPeter Brune if (F) { 14426a388c36SPeter Brune PetscValidPointer(F, 3); 14436a388c36SPeter Brune *F = linesearch->vec_func; 14446a388c36SPeter Brune } 14456a388c36SPeter Brune if (Y) { 14466a388c36SPeter Brune PetscValidPointer(Y, 4); 14476a388c36SPeter Brune *Y = linesearch->vec_update; 14486a388c36SPeter Brune } 14496a388c36SPeter Brune if (W) { 14506a388c36SPeter Brune PetscValidPointer(W, 5); 14516a388c36SPeter Brune *W = linesearch->vec_sol_new; 14526a388c36SPeter Brune } 14536a388c36SPeter Brune if (G) { 14546a388c36SPeter Brune PetscValidPointer(G, 6); 14556a388c36SPeter Brune *G = linesearch->vec_func_new; 14566a388c36SPeter Brune } 14573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14586a388c36SPeter Brune } 14596a388c36SPeter Brune 1460f40b411bSPeter Brune /*@ 1461f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1462f6dfbefdSBarry Smith 1463c3339decSBarry Smith Logically Collective 1464f40b411bSPeter Brune 1465f40b411bSPeter Brune Input Parameters: 146678bcb3b5SPeter Brune + linesearch - linesearch context 14676232e825SPeter Brune . X - Solution vector 14686232e825SPeter Brune . F - Function vector 14696232e825SPeter Brune . Y - Search direction vector 14706232e825SPeter Brune . W - Solution work vector 14716232e825SPeter Brune - G - Function work vector 1472f40b411bSPeter Brune 147378bcb3b5SPeter Brune Level: advanced 1474f40b411bSPeter Brune 1475f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1476f40b411bSPeter Brune @*/ 1477d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1478d71ae5a4SJacob Faibussowitsch { 14796a388c36SPeter Brune PetscFunctionBegin; 1480f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14816a388c36SPeter Brune if (X) { 14826a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 14836a388c36SPeter Brune linesearch->vec_sol = X; 14846a388c36SPeter Brune } 14856a388c36SPeter Brune if (F) { 14866a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 14876a388c36SPeter Brune linesearch->vec_func = F; 14886a388c36SPeter Brune } 14896a388c36SPeter Brune if (Y) { 14906a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 14916a388c36SPeter Brune linesearch->vec_update = Y; 14926a388c36SPeter Brune } 14936a388c36SPeter Brune if (W) { 14946a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 14956a388c36SPeter Brune linesearch->vec_sol_new = W; 14966a388c36SPeter Brune } 14976a388c36SPeter Brune if (G) { 14986a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 14996a388c36SPeter Brune linesearch->vec_func_new = G; 15006a388c36SPeter Brune } 15013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15026a388c36SPeter Brune } 15036a388c36SPeter Brune 1504e7058c64SPeter Brune /*@C 1505f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1506f6dfbefdSBarry Smith `SNESLineSearch` options in the database. 1507e7058c64SPeter Brune 1508c3339decSBarry Smith Logically Collective 1509e7058c64SPeter Brune 1510e7058c64SPeter Brune Input Parameters: 1511f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1512e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1513e7058c64SPeter Brune 1514*20f4b53cSBarry Smith Level: advanced 1515*20f4b53cSBarry Smith 1516f6dfbefdSBarry Smith Note: 1517e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1518e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1519e7058c64SPeter Brune 1520f6dfbefdSBarry Smith .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1521e7058c64SPeter Brune @*/ 1522d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1523d71ae5a4SJacob Faibussowitsch { 1524e7058c64SPeter Brune PetscFunctionBegin; 1525f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15269566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 15273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1528e7058c64SPeter Brune } 1529e7058c64SPeter Brune 1530e7058c64SPeter Brune /*@C 1531f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1532f1c6b773SPeter Brune SNESLineSearch options in the database. 1533e7058c64SPeter Brune 1534e7058c64SPeter Brune Not Collective 1535e7058c64SPeter Brune 1536e7058c64SPeter Brune Input Parameter: 1537f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 1538e7058c64SPeter Brune 1539e7058c64SPeter Brune Output Parameter: 1540e7058c64SPeter Brune . prefix - pointer to the prefix string used 1541e7058c64SPeter Brune 1542e7058c64SPeter Brune Level: advanced 1543e7058c64SPeter Brune 1544*20f4b53cSBarry Smith Fortran Note: 1545*20f4b53cSBarry Smith The user should pass in a string 'prefix' of 1546*20f4b53cSBarry Smith sufficient length to hold the prefix. 1547*20f4b53cSBarry Smith 1548f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1549e7058c64SPeter Brune @*/ 1550d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1551d71ae5a4SJacob Faibussowitsch { 1552e7058c64SPeter Brune PetscFunctionBegin; 1553f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15549566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 15553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1556e7058c64SPeter Brune } 1557bf7f4e0aSPeter Brune 15588d359177SBarry Smith /*@C 1559f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1560f40b411bSPeter Brune 1561d8d19677SJose E. Roman Input Parameters: 1562f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1563f40b411bSPeter Brune - nwork - the number of work vectors 1564f40b411bSPeter Brune 1565f40b411bSPeter Brune Level: developer 1566f40b411bSPeter Brune 1567f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetWorkVecs()` 1568f40b411bSPeter Brune @*/ 1569d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1570d71ae5a4SJacob Faibussowitsch { 1571bf7f4e0aSPeter Brune PetscFunctionBegin; 15720fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 15739566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 15743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1575bf7f4e0aSPeter Brune } 1576bf7f4e0aSPeter Brune 1577f40b411bSPeter Brune /*@ 1578422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1579f40b411bSPeter Brune 1580f40b411bSPeter Brune Input Parameters: 158178bcb3b5SPeter Brune . linesearch - linesearch context 1582f40b411bSPeter Brune 1583f40b411bSPeter Brune Output Parameters: 1584422a814eSBarry Smith . result - The success or failure status 1585f40b411bSPeter Brune 1586*20f4b53cSBarry Smith Level: developer 1587*20f4b53cSBarry Smith 1588f6dfbefdSBarry Smith Note: 1589f6dfbefdSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed 1590cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1591cd7522eaSPeter Brune 1592f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1593f40b411bSPeter Brune @*/ 1594d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1595d71ae5a4SJacob Faibussowitsch { 1596bf7f4e0aSPeter Brune PetscFunctionBegin; 1597f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1598422a814eSBarry Smith PetscValidPointer(result, 2); 1599422a814eSBarry Smith *result = linesearch->result; 16003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1601bf7f4e0aSPeter Brune } 1602bf7f4e0aSPeter Brune 1603f40b411bSPeter Brune /*@ 1604422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1605f40b411bSPeter Brune 1606f40b411bSPeter Brune Input Parameters: 160778bcb3b5SPeter Brune + linesearch - linesearch context 1608422a814eSBarry Smith - result - The success or failure status 1609f40b411bSPeter Brune 1610*20f4b53cSBarry Smith Level: developer 1611*20f4b53cSBarry Smith 1612f6dfbefdSBarry Smith Note: 1613f6dfbefdSBarry Smith This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1614cd7522eaSPeter Brune the success or failure of the line search method. 1615cd7522eaSPeter Brune 1616f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1617f40b411bSPeter Brune @*/ 1618d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1619d71ae5a4SJacob Faibussowitsch { 16206a388c36SPeter Brune PetscFunctionBegin; 16215fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1622422a814eSBarry Smith linesearch->result = result; 16233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16246a388c36SPeter Brune } 16256a388c36SPeter Brune 16269bd66eb0SPeter Brune /*@C 1627f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16289bd66eb0SPeter Brune 1629c3339decSBarry Smith Logically Collective 1630f6dfbefdSBarry Smith 16319bd66eb0SPeter Brune Input Parameters: 1632f6dfbefdSBarry Smith + snes - nonlinear context obtained from `SNESCreate()` 16339bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 16349bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16359bd66eb0SPeter Brune 1636*20f4b53cSBarry Smith Calling sequence of `projectfunc`: 16379bd66eb0SPeter Brune .vb 1638*20f4b53cSBarry Smith PetscErrorCode projectfunc(SNES snes, Vec X) 16399bd66eb0SPeter Brune .ve 16409bd66eb0SPeter Brune + snes - nonlinear context 1641*20f4b53cSBarry Smith - X - current solution, store the projected solution here 16429bd66eb0SPeter Brune 1643*20f4b53cSBarry Smith Calling sequence of `normfunc`: 16449bd66eb0SPeter Brune .vb 1645*20f4b53cSBarry Smith PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm) 16469bd66eb0SPeter Brune .ve 16479bd66eb0SPeter Brune + snes - nonlinear context 16489bd66eb0SPeter Brune . X - current solution 1649*20f4b53cSBarry Smith . F - current residual 1650*20f4b53cSBarry Smith - fnorm - VI-specific norm of the function 16519bd66eb0SPeter Brune 1652f6dfbefdSBarry Smith Level: advanced 16539bd66eb0SPeter Brune 1654*20f4b53cSBarry Smith Note: 1655*20f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 1656*20f4b53cSBarry Smith 1657*20f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated 1658*20f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`. 1659*20f4b53cSBarry Smith 1660f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 16619bd66eb0SPeter Brune @*/ 1662d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1663d71ae5a4SJacob Faibussowitsch { 16649bd66eb0SPeter Brune PetscFunctionBegin; 1665f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16669bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 16679bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 16683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16699bd66eb0SPeter Brune } 16709bd66eb0SPeter Brune 16719bd66eb0SPeter Brune /*@C 1672f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 16739bd66eb0SPeter Brune 1674f6dfbefdSBarry Smith Not Collective 1675f6dfbefdSBarry Smith 1676f899ff85SJose E. Roman Input Parameter: 1677f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()` 16789bd66eb0SPeter Brune 16799bd66eb0SPeter Brune Output Parameters: 16809bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16819bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16829bd66eb0SPeter Brune 1683f6dfbefdSBarry Smith Level: advanced 16849bd66eb0SPeter Brune 1685f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 16869bd66eb0SPeter Brune @*/ 1687d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1688d71ae5a4SJacob Faibussowitsch { 16899bd66eb0SPeter Brune PetscFunctionBegin; 16909bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16919bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16939bd66eb0SPeter Brune } 16949bd66eb0SPeter Brune 1695bf7f4e0aSPeter Brune /*@C 1696f6dfbefdSBarry Smith SNESLineSearchRegister - register a line search method 1697bf7f4e0aSPeter Brune 1698bf7f4e0aSPeter Brune Level: advanced 1699f6dfbefdSBarry Smith 1700f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1701bf7f4e0aSPeter Brune @*/ 1702d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch)) 1703d71ae5a4SJacob Faibussowitsch { 1704bf7f4e0aSPeter Brune PetscFunctionBegin; 17059566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 17069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 17073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1708bf7f4e0aSPeter Brune } 1709