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 142fe279fdSBarry Smith Input Parameter: 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 2220f4b53cSBarry Smith Level: advanced 2320f4b53cSBarry 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 502fe279fdSBarry Smith Input Parameter: 51dcf2fd19SBarry Smith . ls - the linesearch object 52dcf2fd19SBarry Smith 532fe279fdSBarry Smith Level: developer 542fe279fdSBarry Smith 55f6dfbefdSBarry Smith Note: 5620f4b53cSBarry Smith This routine is called by the `SNES` implementations. 57dcf2fd19SBarry Smith It does not typically need to be called by the user. 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 8120f4b53cSBarry Smith monitor routine (use `NULL` if no context is desired) 82dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 8320f4b53cSBarry Smith (may be `NULL`) 8420f4b53cSBarry Smith 8520f4b53cSBarry 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 1512fe279fdSBarry Smith Input Parameter: 152f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context). 153f40b411bSPeter Brune 1542fe279fdSBarry Smith Output Parameter: 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 2102fe279fdSBarry Smith Input Parameter: 211f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 212f40b411bSPeter Brune 21320f4b53cSBarry Smith Level: advanced 21420f4b53cSBarry 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 2452fe279fdSBarry Smith Input Parameter: 246f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 247f40b411bSPeter Brune 24820f4b53cSBarry Smith Level: developer 24920f4b53cSBarry 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 29920f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPreCheck()` 30020f4b53cSBarry 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: 33020f4b53cSBarry Smith + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchPreCheck` 33120f4b53cSBarry 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 35020f4b53cSBarry Smith Logically Collective 35186d74e61SPeter Brune 35286d74e61SPeter Brune Input Parameters: 353f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 35420f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPostCheck` 35520f4b53cSBarry 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` 38520f4b53cSBarry 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 4102fe279fdSBarry Smith Output Parameter: 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 44720f4b53cSBarry Smith Level: developer 44820f4b53cSBarry 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 57120f4b53cSBarry Smith Level: Intermediate 57220f4b53cSBarry 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 6202fe279fdSBarry Smith Input Parameter: 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 65220f4b53cSBarry 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 83420f4b53cSBarry Smith Logically Collective 83520f4b53cSBarry Smith 836f40b411bSPeter Brune Input Parameters: 8372fe279fdSBarry Smith + linesearch - linesearch context 8382fe279fdSBarry Smith - viewer - the viewer to display the line search information 839f40b411bSPeter Brune 840f40b411bSPeter Brune Level: intermediate 841f40b411bSPeter Brune 842f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 843f40b411bSPeter Brune @*/ 844d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 845d71ae5a4SJacob Faibussowitsch { 8467f1410a3SPeter Brune PetscBool iascii; 847bf388a1fSBarry Smith 848bf7f4e0aSPeter Brune PetscFunctionBegin; 8497f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 85048a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 8517f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 8527f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2); 853f40b411bSPeter Brune 8549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 8557f1410a3SPeter Brune if (iascii) { 8569566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 8579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 858dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer); 8599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 8609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 8619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 86263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 8636b2b7091SBarry Smith if (linesearch->ops->precheck) { 8646b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 86563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 8667f1410a3SPeter Brune } else { 86763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 8687f1410a3SPeter Brune } 8697f1410a3SPeter Brune } 87048a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 8717f1410a3SPeter Brune } 8723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 873bf7f4e0aSPeter Brune } 874bf7f4e0aSPeter Brune 875ea5d4fccSPeter Brune /*@C 876a80ff896SJed Brown SNESLineSearchGetType - Gets the linesearch type 877a80ff896SJed Brown 878c3339decSBarry Smith Logically Collective 879a80ff896SJed Brown 8802fe279fdSBarry Smith Input Parameter: 881a80ff896SJed Brown . linesearch - linesearch context 882a80ff896SJed Brown 8832fe279fdSBarry Smith Output Parameter: 8842fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set 885a80ff896SJed Brown 886a80ff896SJed Brown Level: intermediate 887a80ff896SJed Brown 888f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 889a80ff896SJed Brown @*/ 890d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 891d71ae5a4SJacob Faibussowitsch { 892a80ff896SJed Brown PetscFunctionBegin; 893a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 894dadcf809SJacob Faibussowitsch PetscValidPointer(type, 2); 895a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 8963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 897a80ff896SJed Brown } 898a80ff896SJed Brown 899a80ff896SJed Brown /*@C 900f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 901f40b411bSPeter Brune 902c3339decSBarry Smith Logically Collective 903f190f2fcSBarry Smith 904f40b411bSPeter Brune Input Parameters: 9058e557f58SPeter Brune + linesearch - linesearch context 906f40b411bSPeter Brune - type - The type of line search to be used 907f40b411bSPeter Brune 908cd7522eaSPeter Brune Available Types: 909f6dfbefdSBarry Smith + `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step 910f6dfbefdSBarry Smith . `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function 911f6dfbefdSBarry Smith . `SNESLINESEARCHL2` - Secant line search over the L2 norm of the function 912f6dfbefdSBarry Smith . `SNESLINESEARCHCP` - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 913f6dfbefdSBarry Smith . `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch 914f6dfbefdSBarry Smith - `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation 9151fe24845SBarry Smith 9163c7db156SBarry Smith Options Database Key: 9170b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 918cd7522eaSPeter Brune 919f40b411bSPeter Brune Level: intermediate 920f40b411bSPeter Brune 921f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 922f40b411bSPeter Brune @*/ 923d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 924d71ae5a4SJacob Faibussowitsch { 925bf7f4e0aSPeter Brune PetscBool match; 9265f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 927bf7f4e0aSPeter Brune 928bf7f4e0aSPeter Brune PetscFunctionBegin; 929f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 930bf7f4e0aSPeter Brune PetscValidCharPointer(type, 2); 931bf7f4e0aSPeter Brune 9329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 9333ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 934bf7f4e0aSPeter Brune 9359566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 9365f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 937bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 938bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 9399566063dSJacob Faibussowitsch PetscCall((*(linesearch)->ops->destroy)(linesearch)); 9400298fd71SBarry Smith linesearch->ops->destroy = NULL; 941bf7f4e0aSPeter Brune } 942f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9439e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9449e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9459e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9469e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 947bf7f4e0aSPeter Brune 9489566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 9499566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 9503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 951bf7f4e0aSPeter Brune } 952bf7f4e0aSPeter Brune 953f40b411bSPeter Brune /*@ 954f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 955f40b411bSPeter Brune 956f40b411bSPeter Brune Input Parameters: 9578e557f58SPeter Brune + linesearch - linesearch context 958f40b411bSPeter Brune - snes - The snes instance 959f40b411bSPeter Brune 96078bcb3b5SPeter Brune Level: developer 96178bcb3b5SPeter Brune 962f6dfbefdSBarry Smith Note: 963f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 964f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 96578bcb3b5SPeter Brune implementations. 966f40b411bSPeter Brune 967f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 968f40b411bSPeter Brune @*/ 969d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 970d71ae5a4SJacob Faibussowitsch { 971bf7f4e0aSPeter Brune PetscFunctionBegin; 972f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 973bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 974bf7f4e0aSPeter Brune linesearch->snes = snes; 9753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 976bf7f4e0aSPeter Brune } 977bf7f4e0aSPeter Brune 978f40b411bSPeter Brune /*@ 979f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 980f6dfbefdSBarry Smith Having an associated `SNES` is necessary because most line search implementations must be able to 981f6dfbefdSBarry Smith evaluate the function using `SNESComputeFunction()` for the associated `SNES`. This routine 982f6dfbefdSBarry Smith is used in the line search implementations when one must get this associated `SNES` instance. 983f6dfbefdSBarry Smith 984f6dfbefdSBarry Smith Not Collective 985f40b411bSPeter Brune 9862fe279fdSBarry Smith Input Parameter: 9878e557f58SPeter Brune . linesearch - linesearch context 988f40b411bSPeter Brune 9892fe279fdSBarry Smith Output Parameter: 9902fe279fdSBarry Smith . snes - The `SNES` instance 991f40b411bSPeter Brune 9928141a3b9SPeter Brune Level: developer 993f40b411bSPeter Brune 994f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 995f40b411bSPeter Brune @*/ 996d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 997d71ae5a4SJacob Faibussowitsch { 998bf7f4e0aSPeter Brune PetscFunctionBegin; 999f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10006a388c36SPeter Brune PetscValidPointer(snes, 2); 1001bf7f4e0aSPeter Brune *snes = linesearch->snes; 10023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1003bf7f4e0aSPeter Brune } 1004bf7f4e0aSPeter Brune 1005f40b411bSPeter Brune /*@ 1006f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1007f40b411bSPeter Brune 1008f6dfbefdSBarry Smith Not Collective 1009f6dfbefdSBarry Smith 10102fe279fdSBarry Smith Input Parameter: 10118e557f58SPeter Brune . linesearch - linesearch context 1012f40b411bSPeter Brune 10132fe279fdSBarry Smith Output Parameter: 1014f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()` 1015f40b411bSPeter Brune 101678bcb3b5SPeter Brune Level: advanced 101778bcb3b5SPeter Brune 1018f6dfbefdSBarry Smith Note: 10198e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 102078bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 1021f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the 102278bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 102378bcb3b5SPeter Brune 1024f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1025f40b411bSPeter Brune @*/ 1026d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1027d71ae5a4SJacob Faibussowitsch { 10286a388c36SPeter Brune PetscFunctionBegin; 1029f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1030534a8f05SLisandro Dalcin PetscValidRealPointer(lambda, 2); 10316a388c36SPeter Brune *lambda = linesearch->lambda; 10323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10336a388c36SPeter Brune } 10346a388c36SPeter Brune 1035f40b411bSPeter Brune /*@ 1036f6dfbefdSBarry Smith SNESLineSearchSetLambda - Sets the linesearch steplength 1037f40b411bSPeter Brune 1038f40b411bSPeter Brune Input Parameters: 10398e557f58SPeter Brune + linesearch - linesearch context 1040f40b411bSPeter Brune - lambda - The last steplength. 1041f40b411bSPeter Brune 104220f4b53cSBarry Smith Level: advanced 104320f4b53cSBarry Smith 1044f6dfbefdSBarry Smith Note: 1045f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()` 1046f6dfbefdSBarry Smith to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1047cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1048cd7522eaSPeter Brune as an inner scaling parameter. 1049cd7522eaSPeter Brune 1050f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()` 1051f40b411bSPeter Brune @*/ 1052d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1053d71ae5a4SJacob Faibussowitsch { 10546a388c36SPeter Brune PetscFunctionBegin; 1055f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10566a388c36SPeter Brune linesearch->lambda = lambda; 10573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10586a388c36SPeter Brune } 10596a388c36SPeter Brune 1060f40b411bSPeter Brune /*@ 10613c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 106278bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 106378bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 106478bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1065f40b411bSPeter Brune 1066f6dfbefdSBarry Smith Not Collective 1067f6dfbefdSBarry Smith 1068f899ff85SJose E. Roman Input Parameter: 10698e557f58SPeter Brune . linesearch - linesearch context 1070f40b411bSPeter Brune 1071f40b411bSPeter Brune Output Parameters: 1072516fe3c3SPeter Brune + steptol - The minimum steplength 10736cc8e53bSPeter Brune . maxstep - The maximum steplength 1074516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1075516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1076516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1077516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1078f40b411bSPeter Brune 107978bcb3b5SPeter Brune Level: intermediate 108078bcb3b5SPeter Brune 1081f6dfbefdSBarry Smith Note: 108278bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 10833c7d6663SPeter Brune the type requires. 1084516fe3c3SPeter Brune 1085f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1086f40b411bSPeter Brune @*/ 1087d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1088d71ae5a4SJacob Faibussowitsch { 10896a388c36SPeter Brune PetscFunctionBegin; 1090f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1091516fe3c3SPeter Brune if (steptol) { 1092534a8f05SLisandro Dalcin PetscValidRealPointer(steptol, 2); 10936a388c36SPeter Brune *steptol = linesearch->steptol; 1094516fe3c3SPeter Brune } 1095516fe3c3SPeter Brune if (maxstep) { 1096534a8f05SLisandro Dalcin PetscValidRealPointer(maxstep, 3); 1097516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1098516fe3c3SPeter Brune } 1099516fe3c3SPeter Brune if (rtol) { 1100534a8f05SLisandro Dalcin PetscValidRealPointer(rtol, 4); 1101516fe3c3SPeter Brune *rtol = linesearch->rtol; 1102516fe3c3SPeter Brune } 1103516fe3c3SPeter Brune if (atol) { 1104534a8f05SLisandro Dalcin PetscValidRealPointer(atol, 5); 1105516fe3c3SPeter Brune *atol = linesearch->atol; 1106516fe3c3SPeter Brune } 1107516fe3c3SPeter Brune if (ltol) { 1108534a8f05SLisandro Dalcin PetscValidRealPointer(ltol, 6); 1109516fe3c3SPeter Brune *ltol = linesearch->ltol; 1110516fe3c3SPeter Brune } 1111516fe3c3SPeter Brune if (max_its) { 1112534a8f05SLisandro Dalcin PetscValidIntPointer(max_its, 7); 1113516fe3c3SPeter Brune *max_its = linesearch->max_its; 1114516fe3c3SPeter Brune } 11153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11166a388c36SPeter Brune } 11176a388c36SPeter Brune 1118f40b411bSPeter Brune /*@ 11193c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 112078bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 112178bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 112278bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1123f40b411bSPeter Brune 1124f6dfbefdSBarry Smith Collective 1125f6dfbefdSBarry Smith 1126f40b411bSPeter Brune Input Parameters: 11278e557f58SPeter Brune + linesearch - linesearch context 1128516fe3c3SPeter Brune . steptol - The minimum steplength 11296cc8e53bSPeter Brune . maxstep - The maximum steplength 1130516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1131516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1132516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1133516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1134f40b411bSPeter Brune 113520f4b53cSBarry Smith Level: intermediate 113620f4b53cSBarry Smith 1137f6dfbefdSBarry Smith Note: 1138f6dfbefdSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in 113978bcb3b5SPeter Brune place of an argument. 1140f40b411bSPeter Brune 1141f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1142f40b411bSPeter Brune @*/ 1143d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 1144d71ae5a4SJacob Faibussowitsch { 11456a388c36SPeter Brune PetscFunctionBegin; 1146f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1147d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1148d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1149d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1150d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1151d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1152d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch, max_its, 7); 1153d3952184SSatish Balay 115413bcc0bdSJacob Faibussowitsch if (steptol != (PetscReal)PETSC_DEFAULT) { 11555f80ce2aSJacob Faibussowitsch PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 11566a388c36SPeter Brune linesearch->steptol = steptol; 1157d3952184SSatish Balay } 1158d3952184SSatish Balay 115913bcc0bdSJacob Faibussowitsch if (maxstep != (PetscReal)PETSC_DEFAULT) { 11605f80ce2aSJacob Faibussowitsch PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1161516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1162d3952184SSatish Balay } 1163d3952184SSatish Balay 116413bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) { 11652061ca32SJunchao 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); 1166516fe3c3SPeter Brune linesearch->rtol = rtol; 1167d3952184SSatish Balay } 1168d3952184SSatish Balay 116913bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) { 11705f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1171516fe3c3SPeter Brune linesearch->atol = atol; 1172d3952184SSatish Balay } 1173d3952184SSatish Balay 117413bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) { 11755f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1176516fe3c3SPeter Brune linesearch->ltol = ltol; 1177d3952184SSatish Balay } 1178d3952184SSatish Balay 1179d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 118063a3b9bcSJacob Faibussowitsch PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its); 1181516fe3c3SPeter Brune linesearch->max_its = max_its; 1182d3952184SSatish Balay } 11833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11846a388c36SPeter Brune } 11856a388c36SPeter Brune 1186f40b411bSPeter Brune /*@ 1187f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1188f40b411bSPeter Brune 11892fe279fdSBarry Smith Input Parameter: 11908e557f58SPeter Brune . linesearch - linesearch context 1191f40b411bSPeter Brune 11922fe279fdSBarry Smith Output Parameter: 11938e557f58SPeter Brune . damping - The damping parameter 1194f40b411bSPeter Brune 119578bcb3b5SPeter Brune Level: advanced 1196f40b411bSPeter Brune 1197db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN` 1198f40b411bSPeter Brune @*/ 1199f40b411bSPeter Brune 1200d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1201d71ae5a4SJacob Faibussowitsch { 12026a388c36SPeter Brune PetscFunctionBegin; 1203f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1204534a8f05SLisandro Dalcin PetscValidRealPointer(damping, 2); 12056a388c36SPeter Brune *damping = linesearch->damping; 12063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12076a388c36SPeter Brune } 12086a388c36SPeter Brune 1209f40b411bSPeter Brune /*@ 1210fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1211f40b411bSPeter Brune 1212f40b411bSPeter Brune Input Parameters: 121303fd4120SBarry Smith + linesearch - linesearch context 121403fd4120SBarry Smith - damping - The damping parameter 1215f40b411bSPeter Brune 1216f6dfbefdSBarry Smith Options Database Key: 1217f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value 1218f6dfbefdSBarry Smith 1219f40b411bSPeter Brune Level: intermediate 1220f40b411bSPeter Brune 1221f6dfbefdSBarry Smith Note: 1222f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1223cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 122478bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1225cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1226cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1227cd7522eaSPeter Brune 1228f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()` 1229f40b411bSPeter Brune @*/ 1230d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1231d71ae5a4SJacob Faibussowitsch { 12326a388c36SPeter Brune PetscFunctionBegin; 1233f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12346a388c36SPeter Brune linesearch->damping = damping; 12353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12366a388c36SPeter Brune } 12376a388c36SPeter Brune 123859405d5eSPeter Brune /*@ 123959405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 124059405d5eSPeter Brune 1241f6dfbefdSBarry Smith Input Parameter: 124278bcb3b5SPeter Brune . linesearch - linesearch context 124359405d5eSPeter Brune 1244f6dfbefdSBarry Smith Output Parameter: 12458e557f58SPeter Brune . order - The order 124659405d5eSPeter Brune 124778bcb3b5SPeter Brune Possible Values for order: 1248f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1249f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1250f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 125178bcb3b5SPeter Brune 125259405d5eSPeter Brune Level: intermediate 125359405d5eSPeter Brune 1254f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()` 125559405d5eSPeter Brune @*/ 125659405d5eSPeter Brune 1257d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1258d71ae5a4SJacob Faibussowitsch { 125959405d5eSPeter Brune PetscFunctionBegin; 126059405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1261534a8f05SLisandro Dalcin PetscValidIntPointer(order, 2); 126259405d5eSPeter Brune *order = linesearch->order; 12633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126459405d5eSPeter Brune } 126559405d5eSPeter Brune 126659405d5eSPeter Brune /*@ 12671f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 126859405d5eSPeter Brune 126959405d5eSPeter Brune Input Parameters: 1270a2b725a8SWilliam Gropp + linesearch - linesearch context 1271a2b725a8SWilliam Gropp - order - The damping parameter 127259405d5eSPeter Brune 127359405d5eSPeter Brune Level: intermediate 127459405d5eSPeter Brune 127578bcb3b5SPeter Brune Possible Values for order: 1276f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1277f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1278f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 127978bcb3b5SPeter Brune 128059405d5eSPeter Brune Notes: 128159405d5eSPeter Brune Variable orders are supported by the following line searches: 128278bcb3b5SPeter Brune + bt - cubic and quadratic 128378bcb3b5SPeter Brune - cp - linear and quadratic 128459405d5eSPeter Brune 1285f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 128659405d5eSPeter Brune @*/ 1287d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1288d71ae5a4SJacob Faibussowitsch { 128959405d5eSPeter Brune PetscFunctionBegin; 129059405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 129159405d5eSPeter Brune linesearch->order = order; 12923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129359405d5eSPeter Brune } 129459405d5eSPeter Brune 1295f40b411bSPeter Brune /*@ 1296f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1297f40b411bSPeter Brune 1298f6dfbefdSBarry Smith Not Collective 1299f6dfbefdSBarry Smith 1300f899ff85SJose E. Roman Input Parameter: 130178bcb3b5SPeter Brune . linesearch - linesearch context 1302f40b411bSPeter Brune 1303f40b411bSPeter Brune Output Parameters: 1304f40b411bSPeter Brune + xnorm - The norm of the current solution 1305*a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1306f40b411bSPeter Brune - ynorm - The norm of the current update 1307f40b411bSPeter Brune 130878bcb3b5SPeter Brune Level: developer 1309f40b411bSPeter Brune 1310*a68bbae5SBarry Smith Note: 1311*a68bbae5SBarry Smith Some values may not be up to date at particular points in the code. 1312*a68bbae5SBarry Smith 1313*a68bbae5SBarry Smith This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1314*a68bbae5SBarry Smith computed values. 1315*a68bbae5SBarry Smith 1316f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1317f40b411bSPeter Brune @*/ 1318d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1319d71ae5a4SJacob Faibussowitsch { 1320bf7f4e0aSPeter Brune PetscFunctionBegin; 1321f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1322f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1323f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1324f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 13253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1326bf7f4e0aSPeter Brune } 1327bf7f4e0aSPeter Brune 1328f40b411bSPeter Brune /*@ 1329f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1330f40b411bSPeter Brune 1331c3339decSBarry Smith Collective 1332f6dfbefdSBarry Smith 1333f40b411bSPeter Brune Input Parameters: 133478bcb3b5SPeter Brune + linesearch - linesearch context 1335f40b411bSPeter Brune . xnorm - The norm of the current solution 1336*a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1337f40b411bSPeter Brune - ynorm - The norm of the current update 1338f40b411bSPeter Brune 1339f6dfbefdSBarry Smith Level: developer 1340f40b411bSPeter Brune 1341f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1342f40b411bSPeter Brune @*/ 1343d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1344d71ae5a4SJacob Faibussowitsch { 13456a388c36SPeter Brune PetscFunctionBegin; 1346f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 13476a388c36SPeter Brune linesearch->xnorm = xnorm; 13486a388c36SPeter Brune linesearch->fnorm = fnorm; 13496a388c36SPeter Brune linesearch->ynorm = ynorm; 13503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13516a388c36SPeter Brune } 13526a388c36SPeter Brune 1353f40b411bSPeter Brune /*@ 1354f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1355f40b411bSPeter Brune 13562fe279fdSBarry Smith Input Parameter: 135778bcb3b5SPeter Brune . linesearch - linesearch context 1358f40b411bSPeter Brune 135920f4b53cSBarry Smith Options Database Key: 13608e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1361f40b411bSPeter Brune 1362f40b411bSPeter Brune Level: intermediate 1363f40b411bSPeter Brune 1364f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1365f40b411bSPeter Brune @*/ 1366d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1367d71ae5a4SJacob Faibussowitsch { 13689bd66eb0SPeter Brune SNES snes; 1369bf388a1fSBarry Smith 13706a388c36SPeter Brune PetscFunctionBegin; 13716a388c36SPeter Brune if (linesearch->norms) { 13729bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 13739566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 13749566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13759566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13769566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 13779bd66eb0SPeter Brune } else { 13789566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13799566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13809566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13819566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13829566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13839566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13846a388c36SPeter Brune } 13859bd66eb0SPeter Brune } 13863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13876a388c36SPeter Brune } 13886a388c36SPeter Brune 13896f263ca3SPeter Brune /*@ 13906f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 13916f263ca3SPeter Brune 13926f263ca3SPeter Brune Input Parameters: 139378bcb3b5SPeter Brune + linesearch - linesearch context 139478bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 13956f263ca3SPeter Brune 139620f4b53cSBarry Smith Options Database Key: 13970b00b554SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch 13986f263ca3SPeter Brune 139920f4b53cSBarry Smith Level: intermediate 140020f4b53cSBarry Smith 1401f6dfbefdSBarry Smith Note: 1402f6dfbefdSBarry 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. 14036f263ca3SPeter Brune 1404f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 14056f263ca3SPeter Brune @*/ 1406d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1407d71ae5a4SJacob Faibussowitsch { 14086f263ca3SPeter Brune PetscFunctionBegin; 14096f263ca3SPeter Brune linesearch->norms = flg; 14103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14116f263ca3SPeter Brune } 14126f263ca3SPeter Brune 1413f40b411bSPeter Brune /*@ 1414f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1415f6dfbefdSBarry Smith 1416f6dfbefdSBarry Smith Not Collective on but the vectors are parallel 1417f40b411bSPeter Brune 1418f899ff85SJose E. Roman Input Parameter: 141978bcb3b5SPeter Brune . linesearch - linesearch context 1420f40b411bSPeter Brune 1421f40b411bSPeter Brune Output Parameters: 14226232e825SPeter Brune + X - Solution vector 14236232e825SPeter Brune . F - Function vector 14246232e825SPeter Brune . Y - Search direction vector 14256232e825SPeter Brune . W - Solution work vector 14266232e825SPeter Brune - G - Function work vector 14276232e825SPeter Brune 142820f4b53cSBarry Smith Level: advanced 142920f4b53cSBarry Smith 14307bba9028SPeter Brune Notes: 143120f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a 143220f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the 143320f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the 14346232e825SPeter Brune function evaluated at the new solution. 1435f40b411bSPeter Brune 1436f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 14372a7a6963SBarry Smith 1438f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1439f40b411bSPeter Brune @*/ 1440d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1441d71ae5a4SJacob Faibussowitsch { 14426a388c36SPeter Brune PetscFunctionBegin; 1443f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14446a388c36SPeter Brune if (X) { 14456a388c36SPeter Brune PetscValidPointer(X, 2); 14466a388c36SPeter Brune *X = linesearch->vec_sol; 14476a388c36SPeter Brune } 14486a388c36SPeter Brune if (F) { 14496a388c36SPeter Brune PetscValidPointer(F, 3); 14506a388c36SPeter Brune *F = linesearch->vec_func; 14516a388c36SPeter Brune } 14526a388c36SPeter Brune if (Y) { 14536a388c36SPeter Brune PetscValidPointer(Y, 4); 14546a388c36SPeter Brune *Y = linesearch->vec_update; 14556a388c36SPeter Brune } 14566a388c36SPeter Brune if (W) { 14576a388c36SPeter Brune PetscValidPointer(W, 5); 14586a388c36SPeter Brune *W = linesearch->vec_sol_new; 14596a388c36SPeter Brune } 14606a388c36SPeter Brune if (G) { 14616a388c36SPeter Brune PetscValidPointer(G, 6); 14626a388c36SPeter Brune *G = linesearch->vec_func_new; 14636a388c36SPeter Brune } 14643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14656a388c36SPeter Brune } 14666a388c36SPeter Brune 1467f40b411bSPeter Brune /*@ 1468f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1469f6dfbefdSBarry Smith 1470c3339decSBarry Smith Logically Collective 1471f40b411bSPeter Brune 1472f40b411bSPeter Brune Input Parameters: 147378bcb3b5SPeter Brune + linesearch - linesearch context 14746232e825SPeter Brune . X - Solution vector 14756232e825SPeter Brune . F - Function vector 14766232e825SPeter Brune . Y - Search direction vector 14776232e825SPeter Brune . W - Solution work vector 14786232e825SPeter Brune - G - Function work vector 1479f40b411bSPeter Brune 148078bcb3b5SPeter Brune Level: advanced 1481f40b411bSPeter Brune 1482f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1483f40b411bSPeter Brune @*/ 1484d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1485d71ae5a4SJacob Faibussowitsch { 14866a388c36SPeter Brune PetscFunctionBegin; 1487f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14886a388c36SPeter Brune if (X) { 14896a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 14906a388c36SPeter Brune linesearch->vec_sol = X; 14916a388c36SPeter Brune } 14926a388c36SPeter Brune if (F) { 14936a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 14946a388c36SPeter Brune linesearch->vec_func = F; 14956a388c36SPeter Brune } 14966a388c36SPeter Brune if (Y) { 14976a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 14986a388c36SPeter Brune linesearch->vec_update = Y; 14996a388c36SPeter Brune } 15006a388c36SPeter Brune if (W) { 15016a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 15026a388c36SPeter Brune linesearch->vec_sol_new = W; 15036a388c36SPeter Brune } 15046a388c36SPeter Brune if (G) { 15056a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 15066a388c36SPeter Brune linesearch->vec_func_new = G; 15076a388c36SPeter Brune } 15083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15096a388c36SPeter Brune } 15106a388c36SPeter Brune 1511e7058c64SPeter Brune /*@C 1512f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1513f6dfbefdSBarry Smith `SNESLineSearch` options in the database. 1514e7058c64SPeter Brune 1515c3339decSBarry Smith Logically Collective 1516e7058c64SPeter Brune 1517e7058c64SPeter Brune Input Parameters: 1518f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1519e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1520e7058c64SPeter Brune 152120f4b53cSBarry Smith Level: advanced 152220f4b53cSBarry Smith 1523f6dfbefdSBarry Smith Note: 1524e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1525e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1526e7058c64SPeter Brune 1527f6dfbefdSBarry Smith .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1528e7058c64SPeter Brune @*/ 1529d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1530d71ae5a4SJacob Faibussowitsch { 1531e7058c64SPeter Brune PetscFunctionBegin; 1532f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15339566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 15343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1535e7058c64SPeter Brune } 1536e7058c64SPeter Brune 1537e7058c64SPeter Brune /*@C 1538f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1539f1c6b773SPeter Brune SNESLineSearch options in the database. 1540e7058c64SPeter Brune 1541e7058c64SPeter Brune Not Collective 1542e7058c64SPeter Brune 1543e7058c64SPeter Brune Input Parameter: 1544f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 1545e7058c64SPeter Brune 1546e7058c64SPeter Brune Output Parameter: 1547e7058c64SPeter Brune . prefix - pointer to the prefix string used 1548e7058c64SPeter Brune 1549e7058c64SPeter Brune Level: advanced 1550e7058c64SPeter Brune 155120f4b53cSBarry Smith Fortran Note: 155220f4b53cSBarry Smith The user should pass in a string 'prefix' of 155320f4b53cSBarry Smith sufficient length to hold the prefix. 155420f4b53cSBarry Smith 1555f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1556e7058c64SPeter Brune @*/ 1557d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1558d71ae5a4SJacob Faibussowitsch { 1559e7058c64SPeter Brune PetscFunctionBegin; 1560f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 15623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1563e7058c64SPeter Brune } 1564bf7f4e0aSPeter Brune 15658d359177SBarry Smith /*@C 1566f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1567f40b411bSPeter Brune 1568d8d19677SJose E. Roman Input Parameters: 1569f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1570f40b411bSPeter Brune - nwork - the number of work vectors 1571f40b411bSPeter Brune 1572f40b411bSPeter Brune Level: developer 1573f40b411bSPeter Brune 1574f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetWorkVecs()` 1575f40b411bSPeter Brune @*/ 1576d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1577d71ae5a4SJacob Faibussowitsch { 1578bf7f4e0aSPeter Brune PetscFunctionBegin; 15790fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 15809566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 15813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1582bf7f4e0aSPeter Brune } 1583bf7f4e0aSPeter Brune 1584f40b411bSPeter Brune /*@ 1585422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1586f40b411bSPeter Brune 15872fe279fdSBarry Smith Input Parameter: 158878bcb3b5SPeter Brune . linesearch - linesearch context 1589f40b411bSPeter Brune 15902fe279fdSBarry Smith Output Parameter: 1591422a814eSBarry Smith . result - The success or failure status 1592f40b411bSPeter Brune 159320f4b53cSBarry Smith Level: developer 159420f4b53cSBarry Smith 1595f6dfbefdSBarry Smith Note: 1596f6dfbefdSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed 1597cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1598cd7522eaSPeter Brune 1599f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1600f40b411bSPeter Brune @*/ 1601d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1602d71ae5a4SJacob Faibussowitsch { 1603bf7f4e0aSPeter Brune PetscFunctionBegin; 1604f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1605422a814eSBarry Smith PetscValidPointer(result, 2); 1606422a814eSBarry Smith *result = linesearch->result; 16073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1608bf7f4e0aSPeter Brune } 1609bf7f4e0aSPeter Brune 1610f40b411bSPeter Brune /*@ 1611422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1612f40b411bSPeter Brune 1613f40b411bSPeter Brune Input Parameters: 161478bcb3b5SPeter Brune + linesearch - linesearch context 1615422a814eSBarry Smith - result - The success or failure status 1616f40b411bSPeter Brune 161720f4b53cSBarry Smith Level: developer 161820f4b53cSBarry Smith 1619f6dfbefdSBarry Smith Note: 1620f6dfbefdSBarry Smith This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1621cd7522eaSPeter Brune the success or failure of the line search method. 1622cd7522eaSPeter Brune 1623f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1624f40b411bSPeter Brune @*/ 1625d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1626d71ae5a4SJacob Faibussowitsch { 16276a388c36SPeter Brune PetscFunctionBegin; 16285fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1629422a814eSBarry Smith linesearch->result = result; 16303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16316a388c36SPeter Brune } 16326a388c36SPeter Brune 16339bd66eb0SPeter Brune /*@C 1634f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16359bd66eb0SPeter Brune 1636c3339decSBarry Smith Logically Collective 1637f6dfbefdSBarry Smith 16389bd66eb0SPeter Brune Input Parameters: 1639f6dfbefdSBarry Smith + snes - nonlinear context obtained from `SNESCreate()` 16409bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 16419bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16429bd66eb0SPeter Brune 164320f4b53cSBarry Smith Calling sequence of `projectfunc`: 16449bd66eb0SPeter Brune .vb 164520f4b53cSBarry Smith PetscErrorCode projectfunc(SNES snes, Vec X) 16469bd66eb0SPeter Brune .ve 16479bd66eb0SPeter Brune + snes - nonlinear context 164820f4b53cSBarry Smith - X - current solution, store the projected solution here 16499bd66eb0SPeter Brune 165020f4b53cSBarry Smith Calling sequence of `normfunc`: 16519bd66eb0SPeter Brune .vb 165220f4b53cSBarry Smith PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm) 16539bd66eb0SPeter Brune .ve 16549bd66eb0SPeter Brune + snes - nonlinear context 16559bd66eb0SPeter Brune . X - current solution 165620f4b53cSBarry Smith . F - current residual 165720f4b53cSBarry Smith - fnorm - VI-specific norm of the function 16589bd66eb0SPeter Brune 1659f6dfbefdSBarry Smith Level: advanced 16609bd66eb0SPeter Brune 166120f4b53cSBarry Smith Note: 166220f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 166320f4b53cSBarry Smith 166420f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated 166520f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`. 166620f4b53cSBarry Smith 1667f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 16689bd66eb0SPeter Brune @*/ 1669d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1670d71ae5a4SJacob Faibussowitsch { 16719bd66eb0SPeter Brune PetscFunctionBegin; 1672f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16739bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 16749bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 16753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16769bd66eb0SPeter Brune } 16779bd66eb0SPeter Brune 16789bd66eb0SPeter Brune /*@C 1679f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 16809bd66eb0SPeter Brune 1681f6dfbefdSBarry Smith Not Collective 1682f6dfbefdSBarry Smith 1683f899ff85SJose E. Roman Input Parameter: 1684f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()` 16859bd66eb0SPeter Brune 16869bd66eb0SPeter Brune Output Parameters: 16879bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16889bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16899bd66eb0SPeter Brune 1690f6dfbefdSBarry Smith Level: advanced 16919bd66eb0SPeter Brune 1692f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 16939bd66eb0SPeter Brune @*/ 1694d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1695d71ae5a4SJacob Faibussowitsch { 16969bd66eb0SPeter Brune PetscFunctionBegin; 16979bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16989bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17009bd66eb0SPeter Brune } 17019bd66eb0SPeter Brune 1702bf7f4e0aSPeter Brune /*@C 1703f6dfbefdSBarry Smith SNESLineSearchRegister - register a line search method 1704bf7f4e0aSPeter Brune 1705bf7f4e0aSPeter Brune Level: advanced 1706f6dfbefdSBarry Smith 1707f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1708bf7f4e0aSPeter Brune @*/ 1709d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch)) 1710d71ae5a4SJacob Faibussowitsch { 1711bf7f4e0aSPeter Brune PetscFunctionBegin; 17129566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 17139566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 17143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1715bf7f4e0aSPeter Brune } 1716