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 27420bcc1bSBarry Smith This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it 28dcf2fd19SBarry Smith that one. 29dcf2fd19SBarry Smith 30420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `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: 56420bcc1bSBarry Smith This routine is called by the `SNESLineSearch` implementations. 57dcf2fd19SBarry Smith It does not typically need to be called by the user. 58dcf2fd19SBarry Smith 59420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `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 80420bcc1bSBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired) 81420bcc1bSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`) 82420bcc1bSBarry Smith 83420bcc1bSBarry Smith Calling sequence of `f`: 84420bcc1bSBarry Smith + ls - the `SNESLineSearch` context 85420bcc1bSBarry Smith - mctx - [optional] user-defined context for private data for the monitor routine 86420bcc1bSBarry Smith 87420bcc1bSBarry Smith Calling sequence of `monitordestroy`: 88420bcc1bSBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine 8920f4b53cSBarry Smith 9020f4b53cSBarry Smith Level: intermediate 91dcf2fd19SBarry Smith 92f6dfbefdSBarry Smith Note: 93dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 94f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` multiple times; all will be called in the 95dcf2fd19SBarry Smith order in which they were set. 96dcf2fd19SBarry Smith 97420bcc1bSBarry Smith Fortran Note: 98f6dfbefdSBarry Smith Only a single monitor function can be set for each `SNESLineSearch` object 99dcf2fd19SBarry Smith 100420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()` 101dcf2fd19SBarry Smith @*/ 102420bcc1bSBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx)) 103d71ae5a4SJacob Faibussowitsch { 10478064530SBarry Smith PetscInt i; 10578064530SBarry Smith PetscBool identical; 10678064530SBarry Smith 107dcf2fd19SBarry Smith PetscFunctionBegin; 108dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 10978064530SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 1109566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical)); 1113ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 11278064530SBarry Smith } 1135f80ce2aSJacob Faibussowitsch PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 114dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 115dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 116dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void *)mctx; 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118dcf2fd19SBarry Smith } 119dcf2fd19SBarry Smith 120dcf2fd19SBarry Smith /*@C 121f6dfbefdSBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries 122dcf2fd19SBarry Smith 123c3339decSBarry Smith Collective 124dcf2fd19SBarry Smith 125dcf2fd19SBarry Smith Input Parameters: 126420bcc1bSBarry Smith + ls - the `SNESLineSearch` object 127f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat` 128dcf2fd19SBarry Smith 129420bcc1bSBarry Smith Options Database Key: 130420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 131dcf2fd19SBarry Smith 132420bcc1bSBarry Smith Level: developer 133420bcc1bSBarry Smith 134420bcc1bSBarry Smith This is not normally called directly but is passed to `SNESLineSearchMonitorSet()` 135420bcc1bSBarry Smith 136420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()` 137dcf2fd19SBarry Smith @*/ 138d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) 139d71ae5a4SJacob Faibussowitsch { 140d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 141dcf2fd19SBarry Smith Vec Y, W, G; 142dcf2fd19SBarry Smith 143dcf2fd19SBarry Smith PetscFunctionBegin; 1449566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G)); 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n")); 1479566063dSJacob Faibussowitsch PetscCall(VecView(Y, viewer)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n")); 1499566063dSJacob Faibussowitsch PetscCall(VecView(W, viewer)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n")); 1519566063dSJacob Faibussowitsch PetscCall(VecView(G, viewer)); 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 154dcf2fd19SBarry Smith } 155dcf2fd19SBarry Smith 156f40b411bSPeter Brune /*@ 157420bcc1bSBarry Smith SNESLineSearchCreate - Creates a `SNESLineSearch` context. 158f40b411bSPeter Brune 159f6dfbefdSBarry Smith Logically Collective 160f40b411bSPeter Brune 1612fe279fdSBarry Smith Input Parameter: 162f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context). 163f40b411bSPeter Brune 1642fe279fdSBarry Smith Output Parameter: 1658e557f58SPeter Brune . outlinesearch - the new line search context 166f40b411bSPeter Brune 167162e0bf5SPeter Brune Level: developer 168162e0bf5SPeter Brune 169f6dfbefdSBarry Smith Note: 170420bcc1bSBarry Smith The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance 171f6dfbefdSBarry Smith already associated with the `SNES`. 172f40b411bSPeter Brune 173420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()` 174f40b411bSPeter Brune @*/ 175d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 176d71ae5a4SJacob Faibussowitsch { 177f1c6b773SPeter Brune SNESLineSearch linesearch; 178bf388a1fSBarry Smith 179bf7f4e0aSPeter Brune PetscFunctionBegin; 1804f572ea9SToby Isaac PetscAssertPointer(outlinesearch, 2); 1819566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 1820298fd71SBarry Smith *outlinesearch = NULL; 183f5af7f23SKarl Rupp 1849566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView)); 185bf7f4e0aSPeter Brune 1860298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1870298fd71SBarry Smith linesearch->vec_func_new = NULL; 1880298fd71SBarry Smith linesearch->vec_sol = NULL; 1890298fd71SBarry Smith linesearch->vec_func = NULL; 1900298fd71SBarry Smith linesearch->vec_update = NULL; 1919bd66eb0SPeter Brune 192bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 193bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 194bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 195bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 196422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 197bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 198bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 199bf7f4e0aSPeter Brune linesearch->damping = 1.0; 200bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 201bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 202516fe3c3SPeter Brune linesearch->rtol = 1e-8; 203516fe3c3SPeter Brune linesearch->atol = 1e-15; 204516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2050298fd71SBarry Smith linesearch->precheckctx = NULL; 2060298fd71SBarry Smith linesearch->postcheckctx = NULL; 20759405d5eSPeter Brune linesearch->max_its = 1; 208bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2093add74b1SBarry Smith linesearch->monitor = NULL; 210bf7f4e0aSPeter Brune *outlinesearch = linesearch; 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212bf7f4e0aSPeter Brune } 213bf7f4e0aSPeter Brune 214f40b411bSPeter Brune /*@ 21578bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 21678bcb3b5SPeter Brune any required vectors. 217f40b411bSPeter Brune 218c3339decSBarry Smith Collective 219f40b411bSPeter Brune 2202fe279fdSBarry Smith Input Parameter: 221f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 222f40b411bSPeter Brune 22320f4b53cSBarry Smith Level: advanced 22420f4b53cSBarry Smith 225f6dfbefdSBarry Smith Note: 226f6dfbefdSBarry Smith For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`. 227cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 22878bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 229f6dfbefdSBarry Smith of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be 230cd7522eaSPeter Brune allocated upfront. 231cd7522eaSPeter Brune 232420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()` 233f40b411bSPeter Brune @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 235d71ae5a4SJacob Faibussowitsch { 236bf7f4e0aSPeter Brune PetscFunctionBegin; 23748a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 238bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 23948a46eb9SPierre Jolivet if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 24048a46eb9SPierre Jolivet if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 241dbbe0bcdSBarry Smith if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup); 2429566063dSJacob Faibussowitsch if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction)); 243bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 244bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 245bf7f4e0aSPeter Brune } 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247bf7f4e0aSPeter Brune } 248bf7f4e0aSPeter Brune 249f40b411bSPeter Brune /*@ 250f6dfbefdSBarry Smith SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search. 251f40b411bSPeter Brune 252c3339decSBarry Smith Collective 253f40b411bSPeter Brune 2542fe279fdSBarry Smith Input Parameter: 255f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 256f40b411bSPeter Brune 25720f4b53cSBarry Smith Level: developer 25820f4b53cSBarry Smith 259f6dfbefdSBarry Smith Note: 260f6dfbefdSBarry Smith Usually only called by `SNESReset()` 261f190f2fcSBarry Smith 262420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 263f40b411bSPeter Brune @*/ 264d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 265d71ae5a4SJacob Faibussowitsch { 266bf7f4e0aSPeter Brune PetscFunctionBegin; 267dbbe0bcdSBarry Smith if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset); 268f5af7f23SKarl Rupp 2699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_sol_new)); 2709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_func_new)); 271bf7f4e0aSPeter Brune 2729566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 273f5af7f23SKarl Rupp 274bf7f4e0aSPeter Brune linesearch->nwork = 0; 275bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277bf7f4e0aSPeter Brune } 278bf7f4e0aSPeter Brune 279ed07d7d7SPeter Brune /*@C 280f6dfbefdSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search 281f6dfbefdSBarry Smith ` 282e4094ef1SJacob Faibussowitsch 283ed07d7d7SPeter Brune Input Parameters: 284e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context 285e4094ef1SJacob Faibussowitsch - func - function evaluation routine, this is usually the function provided with `SNESSetFunction()` 286ed07d7d7SPeter Brune 287420bcc1bSBarry Smith Calling sequence of `func`: 288420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with 289420bcc1bSBarry Smith . x - the input vector 290420bcc1bSBarry Smith - f - the computed value of the function 291420bcc1bSBarry Smith 292ed07d7d7SPeter Brune Level: developer 293ed07d7d7SPeter Brune 294420bcc1bSBarry Smith Note: 295420bcc1bSBarry Smith By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed 296420bcc1bSBarry Smith 297420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()` 298ed07d7d7SPeter Brune @*/ 299420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f)) 300d71ae5a4SJacob Faibussowitsch { 301ed07d7d7SPeter Brune PetscFunctionBegin; 302ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 303ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 305ed07d7d7SPeter Brune } 306ed07d7d7SPeter Brune 30786d74e61SPeter Brune /*@C 308420bcc1bSBarry Smith SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but 309420bcc1bSBarry Smith before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that 310f190f2fcSBarry Smith determined the search direction. 31186d74e61SPeter Brune 312c3339decSBarry Smith Logically Collective 31386d74e61SPeter Brune 31486d74e61SPeter Brune Input Parameters: 315f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 316420bcc1bSBarry Smith . func - [optional] function evaluation routine 31720f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 31886d74e61SPeter Brune 319420bcc1bSBarry Smith Calling sequence of `func`: 320420bcc1bSBarry Smith + ls - the `SNESLineSearch` context 321420bcc1bSBarry Smith . x - the current solution 322a678f235SPierre Jolivet . d - the current search direction 323420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed 324420bcc1bSBarry Smith - ctx - the context passed to `SNESLineSearchSetPreCheck()` 325420bcc1bSBarry Smith 32686d74e61SPeter Brune Level: intermediate 32786d74e61SPeter Brune 328f6dfbefdSBarry Smith Note: 329420bcc1bSBarry Smith Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete. 330f0b84518SBarry Smith 331f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 332f0b84518SBarry Smith 333420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 334869ba2dcSBarry Smith `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 335f0b84518SBarry Smith 33686d74e61SPeter Brune @*/ 337420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx) 338d71ae5a4SJacob Faibussowitsch { 3399bd66eb0SPeter Brune PetscFunctionBegin; 340f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 341f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 34286d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 3433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34486d74e61SPeter Brune } 34586d74e61SPeter Brune 34686d74e61SPeter Brune /*@C 34753deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 34886d74e61SPeter Brune 349f899ff85SJose E. Roman Input Parameter: 350f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 35186d74e61SPeter Brune 35286d74e61SPeter Brune Output Parameters: 353420bcc1bSBarry Smith + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchSetPreCheck()` 35420f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 35586d74e61SPeter Brune 35686d74e61SPeter Brune Level: intermediate 35786d74e61SPeter Brune 358420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 35986d74e61SPeter Brune @*/ 360d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) 361d71ae5a4SJacob Faibussowitsch { 3629bd66eb0SPeter Brune PetscFunctionBegin; 363f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 364f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 36586d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36786d74e61SPeter Brune } 36886d74e61SPeter Brune 36986d74e61SPeter Brune /*@C 370f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 371f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 37286d74e61SPeter Brune 37320f4b53cSBarry Smith Logically Collective 37486d74e61SPeter Brune 37586d74e61SPeter Brune Input Parameters: 376f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 377420bcc1bSBarry Smith . func - [optional] function evaluation routine 37820f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 37986d74e61SPeter Brune 380420bcc1bSBarry Smith Calling sequence of `func`: 381420bcc1bSBarry Smith + ls - the `SNESLineSearch` context 382420bcc1bSBarry Smith . x - the current solution 383a678f235SPierre Jolivet . d - the current search direction 384a678f235SPierre Jolivet . w - $ w = x + lambda*d $ for some lambda 385420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed 386420bcc1bSBarry Smith . changed_w - indicates `w` has been changed 387420bcc1bSBarry Smith - ctx - the context passed to `SNESLineSearchSetPreCheck()` 388420bcc1bSBarry Smith 38986d74e61SPeter Brune Level: intermediate 39086d74e61SPeter Brune 391f0b84518SBarry Smith Notes: 392420bcc1bSBarry Smith Use `SNESLineSearchSetPreCheck()` to change the step before the line search is complete. 393f0b84518SBarry Smith 394f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 395f0b84518SBarry Smith 396420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`, 397e4094ef1SJacob Faibussowitsch `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 39886d74e61SPeter Brune @*/ 399420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, void *ctx), void *ctx) 400d71ae5a4SJacob Faibussowitsch { 40186d74e61SPeter Brune PetscFunctionBegin; 402f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 403f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 40486d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40686d74e61SPeter Brune } 40786d74e61SPeter Brune 40886d74e61SPeter Brune /*@C 409f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 41086d74e61SPeter Brune 411f899ff85SJose E. Roman Input Parameter: 412f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 41386d74e61SPeter Brune 41486d74e61SPeter Brune Output Parameters: 415420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()` 41620f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 41786d74e61SPeter Brune 41886d74e61SPeter Brune Level: intermediate 41986d74e61SPeter Brune 420420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 42186d74e61SPeter Brune @*/ 422d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 423d71ae5a4SJacob Faibussowitsch { 4249bd66eb0SPeter Brune PetscFunctionBegin; 425f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 426f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 42786d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42986d74e61SPeter Brune } 43086d74e61SPeter Brune 431f40b411bSPeter Brune /*@ 432f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 433f40b411bSPeter Brune 434c3339decSBarry Smith Logically Collective 435f40b411bSPeter Brune 436f40b411bSPeter Brune Input Parameters: 4377b1df9c1SPeter Brune + linesearch - The linesearch instance. 4387b1df9c1SPeter Brune . X - The current solution 4397b1df9c1SPeter Brune - Y - The step direction 440f40b411bSPeter Brune 4412fe279fdSBarry Smith Output Parameter: 442420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y` 443f40b411bSPeter Brune 444f0b84518SBarry Smith Level: advanced 445f40b411bSPeter Brune 446f6dfbefdSBarry Smith Note: 447420bcc1bSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines 448f6dfbefdSBarry Smith 449420bcc1bSBarry Smith Developer Note: 450420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided 451420bcc1bSBarry Smith 452420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, 453fe8e7dddSPierre Jolivet `SNESLineSearchGetPostCheck()` 454f40b411bSPeter Brune @*/ 455d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) 456d71ae5a4SJacob Faibussowitsch { 457bf7f4e0aSPeter Brune PetscFunctionBegin; 458bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 4596b2b7091SBarry Smith if (linesearch->ops->precheck) { 460dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx); 46138bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed, 4); 462bf7f4e0aSPeter Brune } 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 464bf7f4e0aSPeter Brune } 465bf7f4e0aSPeter Brune 466f40b411bSPeter Brune /*@ 467ef46b1a6SStefano Zampini SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 468f40b411bSPeter Brune 469c3339decSBarry Smith Logically Collective 470f40b411bSPeter Brune 471f40b411bSPeter Brune Input Parameters: 4727b1df9c1SPeter Brune + linesearch - The line search context 4737b1df9c1SPeter Brune . X - The last solution 4747b1df9c1SPeter Brune . Y - The step direction 475420bcc1bSBarry Smith - W - The updated solution, $ W = X + lambda*Y $ for some lambda 476f40b411bSPeter Brune 477f40b411bSPeter Brune Output Parameters: 478420bcc1bSBarry Smith + changed_Y - Indicator if the direction `Y` has been changed. 479420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed. 480f40b411bSPeter Brune 48120f4b53cSBarry Smith Level: developer 48220f4b53cSBarry Smith 483f6dfbefdSBarry Smith Note: 484420bcc1bSBarry Smith This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines 485f6dfbefdSBarry Smith 486420bcc1bSBarry Smith Developer Note: 487420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided 488420bcc1bSBarry Smith 489420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 490f40b411bSPeter Brune @*/ 491d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 492d71ae5a4SJacob Faibussowitsch { 493bf7f4e0aSPeter Brune PetscFunctionBegin; 494bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 495bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4966b2b7091SBarry Smith if (linesearch->ops->postcheck) { 497dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx); 49838bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5); 49938bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6); 50086d74e61SPeter Brune } 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50286d74e61SPeter Brune } 50386d74e61SPeter Brune 50486d74e61SPeter Brune /*@C 505*1d27aa22SBarry Smith SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time` 50686d74e61SPeter Brune 507c3339decSBarry Smith Logically Collective 50886d74e61SPeter Brune 5094165533cSJose E. Roman Input Parameters: 510420bcc1bSBarry Smith + linesearch - the line search context 51186d74e61SPeter Brune . X - base state for this step 512907376e6SBarry Smith - ctx - context for this function 51386d74e61SPeter Brune 51497bb3fdcSJose E. Roman Input/Output Parameter: 51597bb3fdcSJose E. Roman . Y - correction, possibly modified 51697bb3fdcSJose E. Roman 51797bb3fdcSJose E. Roman Output Parameter: 518420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified 51986d74e61SPeter Brune 520420bcc1bSBarry Smith Options Database Keys: 521cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 522cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 52386d74e61SPeter Brune 52486d74e61SPeter Brune Level: advanced 52586d74e61SPeter Brune 52686d74e61SPeter Brune Notes: 527f6dfbefdSBarry Smith This function should be passed to `SNESLineSearchSetPreCheck()` 52886d74e61SPeter Brune 52986d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 530*1d27aa22SBarry Smith so the Picard linearization should be provided in place of the "Jacobian" {cite}`hindmarsh1996time`. This correction 53186d74e61SPeter Brune is generally not useful when using a Newton linearization. 53286d74e61SPeter Brune 533420bcc1bSBarry Smith Developer Note: 534420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided 535420bcc1bSBarry Smith 536420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 53786d74e61SPeter Brune @*/ 538d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) 539d71ae5a4SJacob Faibussowitsch { 54086d74e61SPeter Brune PetscReal angle = *(PetscReal *)linesearch->precheckctx; 54186d74e61SPeter Brune Vec Ylast; 54286d74e61SPeter Brune PetscScalar dot; 54386d74e61SPeter Brune PetscInt iter; 54486d74e61SPeter Brune PetscReal ynorm, ylastnorm, theta, angle_radians; 54586d74e61SPeter Brune SNES snes; 54686d74e61SPeter Brune 54786d74e61SPeter Brune PetscFunctionBegin; 5489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 5499566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast)); 55086d74e61SPeter Brune if (!Ylast) { 5519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Ylast)); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Ylast)); 55486d74e61SPeter Brune } 5559566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 55686d74e61SPeter Brune if (iter < 2) { 5579566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 55886d74e61SPeter Brune *changed = PETSC_FALSE; 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56086d74e61SPeter Brune } 56186d74e61SPeter Brune 5629566063dSJacob Faibussowitsch PetscCall(VecDot(Y, Ylast, &dot)); 5639566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); 5649566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 56586d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 566255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0)); 56786d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 56886d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 56986d74e61SPeter Brune /* Modify the step Y */ 57086d74e61SPeter Brune PetscReal alpha, ydiffnorm; 5719566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ylast, -1.0, Y)); 5729566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm)); 573f85e2ce2SBarry Smith alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0; 5749566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 5759566063dSJacob Faibussowitsch PetscCall(VecScale(Y, alpha)); 5769566063dSJacob 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)); 577a47ec85fSBarry Smith *changed = PETSC_TRUE; 57886d74e61SPeter Brune } else { 5799566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle)); 5809566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 58186d74e61SPeter Brune *changed = PETSC_FALSE; 582bf7f4e0aSPeter Brune } 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 584bf7f4e0aSPeter Brune } 585bf7f4e0aSPeter Brune 586f40b411bSPeter Brune /*@ 587cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 588f40b411bSPeter Brune 589c3339decSBarry Smith Collective 590f40b411bSPeter Brune 591f40b411bSPeter Brune Input Parameters: 5928e557f58SPeter Brune + linesearch - The line search context 5938e557f58SPeter Brune - Y - The search direction 594f40b411bSPeter Brune 5956b867d5aSJose E. Roman Input/Output Parameters: 5966b867d5aSJose E. Roman + X - The current solution, on output the new solution 597420bcc1bSBarry Smith . F - The current function value, on output the new function value at the solution value `X` 598420bcc1bSBarry Smith - fnorm - The current norm of `F`, on output the new norm of `F` 599f40b411bSPeter Brune 600cd7522eaSPeter Brune Options Database Keys: 6010b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell 602dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 6031fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 6041fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 6053c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 6063c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 607cd7522eaSPeter Brune 608e4094ef1SJacob Faibussowitsch Level: intermediate 60920f4b53cSBarry Smith 610cd7522eaSPeter Brune Notes: 611f6dfbefdSBarry Smith This is typically called from within a `SNESSolve()` implementation in order to 612f6dfbefdSBarry Smith help with convergence of the nonlinear method. Various `SNES` types use line searches 613cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 614cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 615f6dfbefdSBarry Smith application of the line search may invoke `SNESComputeFunction()` several times, and 616cd7522eaSPeter Brune therefore may be fairly expensive. 617cd7522eaSPeter Brune 618420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`, 619db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetType()` 620f40b411bSPeter Brune @*/ 621d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) 622d71ae5a4SJacob Faibussowitsch { 623bf388a1fSBarry Smith PetscFunctionBegin; 624f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 625bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 626bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 627064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Y, VEC_CLASSID, 5); 628bf7f4e0aSPeter Brune 629422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 630bf7f4e0aSPeter Brune 631bf7f4e0aSPeter Brune linesearch->vec_sol = X; 632bf7f4e0aSPeter Brune linesearch->vec_update = Y; 633bf7f4e0aSPeter Brune linesearch->vec_func = F; 634bf7f4e0aSPeter Brune 6359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(linesearch)); 636bf7f4e0aSPeter Brune 637f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 638bf7f4e0aSPeter Brune 639f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 64048a46eb9SPierre Jolivet else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 641bf7f4e0aSPeter Brune 6429566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 643bf7f4e0aSPeter Brune 644dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, apply); 645bf7f4e0aSPeter Brune 6469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 647bf7f4e0aSPeter Brune 648f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 6493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 650bf7f4e0aSPeter Brune } 651bf7f4e0aSPeter Brune 652f40b411bSPeter Brune /*@ 653f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 654f40b411bSPeter Brune 655c3339decSBarry Smith Collective 656f40b411bSPeter Brune 6572fe279fdSBarry Smith Input Parameter: 6588e557f58SPeter Brune . linesearch - The line search context 659f40b411bSPeter Brune 66084238204SBarry Smith Level: developer 661f40b411bSPeter Brune 662420bcc1bSBarry Smith Note: 663420bcc1bSBarry Smith The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed 664420bcc1bSBarry Smith 665420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()` 666f40b411bSPeter Brune @*/ 667d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) 668d71ae5a4SJacob Faibussowitsch { 669bf7f4e0aSPeter Brune PetscFunctionBegin; 6703ba16761SJacob Faibussowitsch if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS); 671f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1); 6729371c9d4SSatish Balay if (--((PetscObject)(*linesearch))->refct > 0) { 6739371c9d4SSatish Balay *linesearch = NULL; 6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6759371c9d4SSatish Balay } 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch)); 6779566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset(*linesearch)); 6783ba16761SJacob Faibussowitsch PetscTryTypeMethod(*linesearch, destroy); 6799566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 6809566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorCancel((*linesearch))); 6819566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(linesearch)); 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 683bf7f4e0aSPeter Brune } 684bf7f4e0aSPeter Brune 685f40b411bSPeter Brune /*@ 686dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 687bf7f4e0aSPeter Brune 688c3339decSBarry Smith Logically Collective 689f6dfbefdSBarry Smith 690bf7f4e0aSPeter Brune Input Parameters: 691dcf2fd19SBarry Smith + linesearch - the linesearch object 69220f4b53cSBarry Smith - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor 693bf7f4e0aSPeter Brune 694f6dfbefdSBarry Smith Options Database Key: 695dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 696bf7f4e0aSPeter Brune 697bf7f4e0aSPeter Brune Level: intermediate 698bf7f4e0aSPeter Brune 699e4094ef1SJacob Faibussowitsch Developer Notes: 700f6dfbefdSBarry Smith This monitor is implemented differently than the other line search monitors that are set with 701f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the 702d12e167eSBarry Smith line search that are not visible to the other monitors. 703dcf2fd19SBarry Smith 704420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`, 705f6dfbefdSBarry Smith `SNESLineSearchMonitorSetFromOptions()` 706bf7f4e0aSPeter Brune @*/ 707d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 708d71ae5a4SJacob Faibussowitsch { 709bf7f4e0aSPeter Brune PetscFunctionBegin; 7109566063dSJacob Faibussowitsch if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer)); 7119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&linesearch->monitor)); 712dcf2fd19SBarry Smith linesearch->monitor = viewer; 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 714bf7f4e0aSPeter Brune } 715bf7f4e0aSPeter Brune 716f40b411bSPeter Brune /*@ 717420bcc1bSBarry Smith SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()` 718f6dfbefdSBarry Smith 719c3339decSBarry Smith Logically Collective 7206a388c36SPeter Brune 721f190f2fcSBarry Smith Input Parameter: 722420bcc1bSBarry Smith . linesearch - the line search context 723f40b411bSPeter Brune 724f190f2fcSBarry Smith Output Parameter: 7258e557f58SPeter Brune . monitor - monitor context 726f40b411bSPeter Brune 727f40b411bSPeter Brune Level: intermediate 728f40b411bSPeter Brune 729420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 730f40b411bSPeter Brune @*/ 731d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 732d71ae5a4SJacob Faibussowitsch { 7336a388c36SPeter Brune PetscFunctionBegin; 734f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 7356a388c36SPeter Brune *monitor = linesearch->monitor; 7363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7376a388c36SPeter Brune } 7386a388c36SPeter Brune 739dcf2fd19SBarry Smith /*@C 740420bcc1bSBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database 741dcf2fd19SBarry Smith 742c3339decSBarry Smith Collective 743dcf2fd19SBarry Smith 744dcf2fd19SBarry Smith Input Parameters: 745420bcc1bSBarry Smith + ls - `SNESLineSearch` object to monitor 746420bcc1bSBarry Smith . name - the monitor type 747dcf2fd19SBarry Smith . help - message indicating what monitoring is done 748dcf2fd19SBarry Smith . manual - manual page for the monitor 749dcf2fd19SBarry Smith . monitor - the monitor function 750f6dfbefdSBarry 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` 751dcf2fd19SBarry Smith 752420bcc1bSBarry Smith Calling sequence of `monitor`: 753420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored 754420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used 755dcf2fd19SBarry Smith 756420bcc1bSBarry Smith Calling sequence of `monitorsetup`: 757420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored 758420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used 759420bcc1bSBarry Smith 760420bcc1bSBarry Smith Level: advanced 761420bcc1bSBarry Smith 762420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 763db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 764e4094ef1SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 765db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 766c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 767db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 768db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 769dcf2fd19SBarry Smith @*/ 770420bcc1bSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch ls, PetscViewerAndFormat *vf), PetscErrorCode (*monitorsetup)(SNESLineSearch ls, PetscViewerAndFormat *vf)) 771d71ae5a4SJacob Faibussowitsch { 772dcf2fd19SBarry Smith PetscViewer viewer; 773dcf2fd19SBarry Smith PetscViewerFormat format; 774dcf2fd19SBarry Smith PetscBool flg; 775dcf2fd19SBarry Smith 776dcf2fd19SBarry Smith PetscFunctionBegin; 7779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg)); 778dcf2fd19SBarry Smith if (flg) { 779d12e167eSBarry Smith PetscViewerAndFormat *vf; 7809566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 7819566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 7821baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls, vf)); 7839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 784dcf2fd19SBarry Smith } 7853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 786dcf2fd19SBarry Smith } 787dcf2fd19SBarry Smith 788f40b411bSPeter Brune /*@ 789f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 790f40b411bSPeter Brune 791c3339decSBarry Smith Logically Collective 792f6dfbefdSBarry Smith 793f899ff85SJose E. Roman Input Parameter: 794420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context 795f40b411bSPeter Brune 796f40b411bSPeter Brune Options Database Keys: 7970b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 7983c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 799f6dfbefdSBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`) 80071eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 8011a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 8021a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 8031a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 8041a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 8051a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 806dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 807dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 8088e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 809cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 8101a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 811d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 812f40b411bSPeter Brune 813f40b411bSPeter Brune Level: intermediate 814f40b411bSPeter Brune 815420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 816db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 817f40b411bSPeter Brune @*/ 818d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 819d71ae5a4SJacob Faibussowitsch { 8201a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 821bf7f4e0aSPeter Brune char type[256]; 822bf7f4e0aSPeter Brune PetscBool flg, set; 823dcf2fd19SBarry Smith PetscViewer viewer; 824bf388a1fSBarry Smith 825bf7f4e0aSPeter Brune PetscFunctionBegin; 8269566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll()); 827bf7f4e0aSPeter Brune 828d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch); 829f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 8309566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg)); 831bf7f4e0aSPeter Brune if (flg) { 8329566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, type)); 833bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 8349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, deft)); 835bf7f4e0aSPeter Brune } 836bf7f4e0aSPeter Brune 8379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set)); 838dcf2fd19SBarry Smith if (set) { 8399566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer)); 8409566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 841dcf2fd19SBarry Smith } 8429566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL)); 843bf7f4e0aSPeter Brune 8441a9b3a06SPeter Brune /* tolerances */ 8459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL)); 8469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL)); 8479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL)); 8489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL)); 8499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL)); 8509566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL)); 8517a35526eSPeter Brune 8521a9b3a06SPeter Brune /* damping parameters */ 8539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 8541a9b3a06SPeter Brune 8559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 8561a9b3a06SPeter Brune 8571a9b3a06SPeter Brune /* precheck */ 8589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 8597a35526eSPeter Brune if (set) { 8607a35526eSPeter Brune if (flg) { 8617a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 862f5af7f23SKarl Rupp 863d0609cedSBarry 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)); 8649566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 8657a35526eSPeter Brune } else { 8669566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 8677a35526eSPeter Brune } 8687a35526eSPeter Brune } 8699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 8709566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 8717a35526eSPeter Brune 872dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 8735a9b6599SPeter Brune 874dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 875d0609cedSBarry Smith PetscOptionsEnd(); 8763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 877bf7f4e0aSPeter Brune } 878bf7f4e0aSPeter Brune 879f40b411bSPeter Brune /*@ 880f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 881f40b411bSPeter Brune 88220f4b53cSBarry Smith Logically Collective 88320f4b53cSBarry Smith 884f40b411bSPeter Brune Input Parameters: 8852fe279fdSBarry Smith + linesearch - line search context 886420bcc1bSBarry Smith - viewer - the `PetscViewer` to display the line search information to 887f40b411bSPeter Brune 888f40b411bSPeter Brune Level: intermediate 889f40b411bSPeter Brune 890420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 891f40b411bSPeter Brune @*/ 892d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 893d71ae5a4SJacob Faibussowitsch { 8947f1410a3SPeter Brune PetscBool iascii; 895bf388a1fSBarry Smith 896bf7f4e0aSPeter Brune PetscFunctionBegin; 8977f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 89848a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 8997f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 9007f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2); 901f40b411bSPeter Brune 9029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 9037f1410a3SPeter Brune if (iascii) { 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 9059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 906dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer); 9079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 9089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 9099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 91063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 9116b2b7091SBarry Smith if (linesearch->ops->precheck) { 9126b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 91363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 9147f1410a3SPeter Brune } else { 91563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 9167f1410a3SPeter Brune } 9177f1410a3SPeter Brune } 91848a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 9197f1410a3SPeter Brune } 9203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 921bf7f4e0aSPeter Brune } 922bf7f4e0aSPeter Brune 923ea5d4fccSPeter Brune /*@C 924420bcc1bSBarry Smith SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch` 925a80ff896SJed Brown 926c3339decSBarry Smith Logically Collective 927a80ff896SJed Brown 9282fe279fdSBarry Smith Input Parameter: 929420bcc1bSBarry Smith . linesearch - the line search context 930a80ff896SJed Brown 9312fe279fdSBarry Smith Output Parameter: 9322fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set 933a80ff896SJed Brown 934a80ff896SJed Brown Level: intermediate 935a80ff896SJed Brown 936420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 937a80ff896SJed Brown @*/ 938d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 939d71ae5a4SJacob Faibussowitsch { 940a80ff896SJed Brown PetscFunctionBegin; 941a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9424f572ea9SToby Isaac PetscAssertPointer(type, 2); 943a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 9443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 945a80ff896SJed Brown } 946a80ff896SJed Brown 947a80ff896SJed Brown /*@C 948420bcc1bSBarry Smith SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` 949f40b411bSPeter Brune 950c3339decSBarry Smith Logically Collective 951f190f2fcSBarry Smith 952f40b411bSPeter Brune Input Parameters: 953420bcc1bSBarry Smith + linesearch - the line search context 954ceaaa498SBarry Smith - type - The type of line search to be used, see `SNESLineSearchType` 9551fe24845SBarry Smith 9563c7db156SBarry Smith Options Database Key: 9570b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 958cd7522eaSPeter Brune 959f40b411bSPeter Brune Level: intermediate 960f40b411bSPeter Brune 961420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 962f40b411bSPeter Brune @*/ 963d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 964d71ae5a4SJacob Faibussowitsch { 965bf7f4e0aSPeter Brune PetscBool match; 9665f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 967bf7f4e0aSPeter Brune 968bf7f4e0aSPeter Brune PetscFunctionBegin; 969f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9704f572ea9SToby Isaac PetscAssertPointer(type, 2); 971bf7f4e0aSPeter Brune 9729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 9733ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 974bf7f4e0aSPeter Brune 9759566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 9766adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 977bf7f4e0aSPeter Brune /* Destroy the previous private line search context */ 978bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 9799566063dSJacob Faibussowitsch PetscCall((*(linesearch)->ops->destroy)(linesearch)); 9800298fd71SBarry Smith linesearch->ops->destroy = NULL; 981bf7f4e0aSPeter Brune } 982f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9839e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9849e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9859e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9869e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 987bf7f4e0aSPeter Brune 9889566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 9899566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 991bf7f4e0aSPeter Brune } 992bf7f4e0aSPeter Brune 993f40b411bSPeter Brune /*@ 994f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 995f40b411bSPeter Brune 996f40b411bSPeter Brune Input Parameters: 997420bcc1bSBarry Smith + linesearch - the line search context 998420bcc1bSBarry Smith - snes - The `SNES` instance 999f40b411bSPeter Brune 100078bcb3b5SPeter Brune Level: developer 100178bcb3b5SPeter Brune 1002f6dfbefdSBarry Smith Note: 1003f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 1004f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 100578bcb3b5SPeter Brune implementations. 1006f40b411bSPeter Brune 1007420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()` 1008f40b411bSPeter Brune @*/ 1009d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1010d71ae5a4SJacob Faibussowitsch { 1011bf7f4e0aSPeter Brune PetscFunctionBegin; 1012f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1013bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 1014bf7f4e0aSPeter Brune linesearch->snes = snes; 10153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1016bf7f4e0aSPeter Brune } 1017bf7f4e0aSPeter Brune 1018f40b411bSPeter Brune /*@ 1019f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 1020f6dfbefdSBarry Smith 1021f6dfbefdSBarry Smith Not Collective 1022f40b411bSPeter Brune 10232fe279fdSBarry Smith Input Parameter: 1024420bcc1bSBarry Smith . linesearch - the line search context 1025f40b411bSPeter Brune 10262fe279fdSBarry Smith Output Parameter: 10272fe279fdSBarry Smith . snes - The `SNES` instance 1028f40b411bSPeter Brune 10298141a3b9SPeter Brune Level: developer 1030f40b411bSPeter Brune 1031420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()` 1032f40b411bSPeter Brune @*/ 1033d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1034d71ae5a4SJacob Faibussowitsch { 1035bf7f4e0aSPeter Brune PetscFunctionBegin; 1036f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10374f572ea9SToby Isaac PetscAssertPointer(snes, 2); 1038bf7f4e0aSPeter Brune *snes = linesearch->snes; 10393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1040bf7f4e0aSPeter Brune } 1041bf7f4e0aSPeter Brune 1042f40b411bSPeter Brune /*@ 1043420bcc1bSBarry Smith SNESLineSearchGetLambda - Gets the last line search steplength used 1044f40b411bSPeter Brune 1045f6dfbefdSBarry Smith Not Collective 1046f6dfbefdSBarry Smith 10472fe279fdSBarry Smith Input Parameter: 1048420bcc1bSBarry Smith . linesearch - the line search context 1049f40b411bSPeter Brune 10502fe279fdSBarry Smith Output Parameter: 1051f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()` 1052f40b411bSPeter Brune 105378bcb3b5SPeter Brune Level: advanced 105478bcb3b5SPeter Brune 1055f6dfbefdSBarry Smith Note: 10568e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 105778bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 1058f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the 105978bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 106078bcb3b5SPeter Brune 1061420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1062f40b411bSPeter Brune @*/ 1063d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1064d71ae5a4SJacob Faibussowitsch { 10656a388c36SPeter Brune PetscFunctionBegin; 1066f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10674f572ea9SToby Isaac PetscAssertPointer(lambda, 2); 10686a388c36SPeter Brune *lambda = linesearch->lambda; 10693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10706a388c36SPeter Brune } 10716a388c36SPeter Brune 1072f40b411bSPeter Brune /*@ 1073f6dfbefdSBarry Smith SNESLineSearchSetLambda - Sets the line search steplength 1074f40b411bSPeter Brune 1075f40b411bSPeter Brune Input Parameters: 10768e557f58SPeter Brune + linesearch - line search context 1077420bcc1bSBarry Smith - lambda - The steplength to use 1078f40b411bSPeter Brune 107920f4b53cSBarry Smith Level: advanced 108020f4b53cSBarry Smith 1081f6dfbefdSBarry Smith Note: 1082f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()` 1083f6dfbefdSBarry Smith to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1084cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1085cd7522eaSPeter Brune as an inner scaling parameter. 1086cd7522eaSPeter Brune 1087420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()` 1088f40b411bSPeter Brune @*/ 1089d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1090d71ae5a4SJacob Faibussowitsch { 10916a388c36SPeter Brune PetscFunctionBegin; 1092f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10936a388c36SPeter Brune linesearch->lambda = lambda; 10943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10956a388c36SPeter Brune } 10966a388c36SPeter Brune 1097f40b411bSPeter Brune /*@ 1098ceaaa498SBarry Smith SNESLineSearchGetTolerances - Gets the tolerances for the line search. 1099f40b411bSPeter Brune 1100f6dfbefdSBarry Smith Not Collective 1101f6dfbefdSBarry Smith 1102f899ff85SJose E. Roman Input Parameter: 1103420bcc1bSBarry Smith . linesearch - the line search context 1104f40b411bSPeter Brune 1105f40b411bSPeter Brune Output Parameters: 1106b13b64c2SBarry Smith + steptol - The minimum steplength 11076cc8e53bSPeter Brune . maxstep - The maximum steplength 1108516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1109516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1110516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1111e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search 1112f40b411bSPeter Brune 111378bcb3b5SPeter Brune Level: intermediate 111478bcb3b5SPeter Brune 1115f6dfbefdSBarry Smith Note: 111678bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 11173c7d6663SPeter Brune the type requires. 1118516fe3c3SPeter Brune 1119420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1120f40b411bSPeter Brune @*/ 1121d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1122d71ae5a4SJacob Faibussowitsch { 11236a388c36SPeter Brune PetscFunctionBegin; 1124f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1125516fe3c3SPeter Brune if (steptol) { 11264f572ea9SToby Isaac PetscAssertPointer(steptol, 2); 11276a388c36SPeter Brune *steptol = linesearch->steptol; 1128516fe3c3SPeter Brune } 1129516fe3c3SPeter Brune if (maxstep) { 11304f572ea9SToby Isaac PetscAssertPointer(maxstep, 3); 1131b13b64c2SBarry Smith *maxstep = linesearch->maxstep; 1132516fe3c3SPeter Brune } 1133516fe3c3SPeter Brune if (rtol) { 11344f572ea9SToby Isaac PetscAssertPointer(rtol, 4); 1135516fe3c3SPeter Brune *rtol = linesearch->rtol; 1136516fe3c3SPeter Brune } 1137516fe3c3SPeter Brune if (atol) { 11384f572ea9SToby Isaac PetscAssertPointer(atol, 5); 1139516fe3c3SPeter Brune *atol = linesearch->atol; 1140516fe3c3SPeter Brune } 1141516fe3c3SPeter Brune if (ltol) { 11424f572ea9SToby Isaac PetscAssertPointer(ltol, 6); 1143516fe3c3SPeter Brune *ltol = linesearch->ltol; 1144516fe3c3SPeter Brune } 1145516fe3c3SPeter Brune if (max_its) { 11464f572ea9SToby Isaac PetscAssertPointer(max_its, 7); 1147516fe3c3SPeter Brune *max_its = linesearch->max_its; 1148516fe3c3SPeter Brune } 11493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11506a388c36SPeter Brune } 11516a388c36SPeter Brune 1152f40b411bSPeter Brune /*@ 1153ceaaa498SBarry Smith SNESLineSearchSetTolerances - Sets the tolerances for the linesearch. 1154f40b411bSPeter Brune 1155f6dfbefdSBarry Smith Collective 1156f6dfbefdSBarry Smith 1157f40b411bSPeter Brune Input Parameters: 1158420bcc1bSBarry Smith + linesearch - the line search context 1159516fe3c3SPeter Brune . steptol - The minimum steplength 11606cc8e53bSPeter Brune . maxstep - The maximum steplength 1161516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1162516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1163516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1164420bcc1bSBarry Smith - max_it - The maximum number of iterations of the line search 1165420bcc1bSBarry Smith 1166420bcc1bSBarry Smith Options Database Keys: 1167420bcc1bSBarry Smith + -snes_linesearch_minlambda - The minimum step length 1168420bcc1bSBarry Smith . -snes_linesearch_maxstep - The maximum step size 1169420bcc1bSBarry Smith . -snes_linesearch_rtol - Relative tolerance for iterative line searches 1170420bcc1bSBarry Smith . -snes_linesearch_atol - Absolute tolerance for iterative line searches 1171420bcc1bSBarry Smith . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 1172420bcc1bSBarry Smith - -snes_linesearch_max_it - The number of iterations for iterative line searches 1173f40b411bSPeter Brune 117420f4b53cSBarry Smith Level: intermediate 117520f4b53cSBarry Smith 1176f6dfbefdSBarry Smith Note: 1177420bcc1bSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument. 1178f40b411bSPeter Brune 1179420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1180f40b411bSPeter Brune @*/ 1181420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it) 1182d71ae5a4SJacob Faibussowitsch { 11836a388c36SPeter Brune PetscFunctionBegin; 1184f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1185d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1186d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1187d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1188d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1189d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1190420bcc1bSBarry Smith PetscValidLogicalCollectiveInt(linesearch, max_it, 7); 1191d3952184SSatish Balay 119213bcc0bdSJacob Faibussowitsch if (steptol != (PetscReal)PETSC_DEFAULT) { 11935f80ce2aSJacob Faibussowitsch PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 11946a388c36SPeter Brune linesearch->steptol = steptol; 1195d3952184SSatish Balay } 1196d3952184SSatish Balay 119713bcc0bdSJacob Faibussowitsch if (maxstep != (PetscReal)PETSC_DEFAULT) { 11985f80ce2aSJacob Faibussowitsch PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1199516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1200d3952184SSatish Balay } 1201d3952184SSatish Balay 120213bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) { 12032061ca32SJunchao 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); 1204516fe3c3SPeter Brune linesearch->rtol = rtol; 1205d3952184SSatish Balay } 1206d3952184SSatish Balay 120713bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) { 12085f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1209516fe3c3SPeter Brune linesearch->atol = atol; 1210d3952184SSatish Balay } 1211d3952184SSatish Balay 121213bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) { 12135f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1214516fe3c3SPeter Brune linesearch->ltol = ltol; 1215d3952184SSatish Balay } 1216d3952184SSatish Balay 1217420bcc1bSBarry Smith if (max_it != PETSC_DEFAULT) { 1218420bcc1bSBarry Smith PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it); 1219420bcc1bSBarry Smith linesearch->max_its = max_it; 1220d3952184SSatish Balay } 12213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12226a388c36SPeter Brune } 12236a388c36SPeter Brune 1224f40b411bSPeter Brune /*@ 1225f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1226f40b411bSPeter Brune 12272fe279fdSBarry Smith Input Parameter: 1228420bcc1bSBarry Smith . linesearch - the line search context 1229f40b411bSPeter Brune 12302fe279fdSBarry Smith Output Parameter: 12318e557f58SPeter Brune . damping - The damping parameter 1232f40b411bSPeter Brune 123378bcb3b5SPeter Brune Level: advanced 1234f40b411bSPeter Brune 1235420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN` 1236f40b411bSPeter Brune @*/ 1237d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1238d71ae5a4SJacob Faibussowitsch { 12396a388c36SPeter Brune PetscFunctionBegin; 1240f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12414f572ea9SToby Isaac PetscAssertPointer(damping, 2); 12426a388c36SPeter Brune *damping = linesearch->damping; 12433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12446a388c36SPeter Brune } 12456a388c36SPeter Brune 1246f40b411bSPeter Brune /*@ 1247fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1248f40b411bSPeter Brune 1249f40b411bSPeter Brune Input Parameters: 1250420bcc1bSBarry Smith + linesearch - the line search context 125103fd4120SBarry Smith - damping - The damping parameter 1252f40b411bSPeter Brune 1253f6dfbefdSBarry Smith Options Database Key: 1254f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value 1255f6dfbefdSBarry Smith 1256f40b411bSPeter Brune Level: intermediate 1257f40b411bSPeter Brune 1258f6dfbefdSBarry Smith Note: 1259f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1260420bcc1bSBarry Smith The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle; 126178bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1262420bcc1bSBarry Smith step may be of greater length than the damping parameter. In the `SNESLINESEARCHBT` line search it is 1263420bcc1bSBarry Smith used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks. 1264cd7522eaSPeter Brune 1265420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()` 1266f40b411bSPeter Brune @*/ 1267d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1268d71ae5a4SJacob Faibussowitsch { 12696a388c36SPeter Brune PetscFunctionBegin; 1270f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12716a388c36SPeter Brune linesearch->damping = damping; 12723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12736a388c36SPeter Brune } 12746a388c36SPeter Brune 127559405d5eSPeter Brune /*@ 127659405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 127759405d5eSPeter Brune 1278f6dfbefdSBarry Smith Input Parameter: 1279420bcc1bSBarry Smith . linesearch - the line search context 128059405d5eSPeter Brune 1281f6dfbefdSBarry Smith Output Parameter: 12828e557f58SPeter Brune . order - The order 128359405d5eSPeter Brune 128459405d5eSPeter Brune Level: intermediate 128559405d5eSPeter Brune 1286420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()` 128759405d5eSPeter Brune @*/ 1288d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1289d71ae5a4SJacob Faibussowitsch { 129059405d5eSPeter Brune PetscFunctionBegin; 129159405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12924f572ea9SToby Isaac PetscAssertPointer(order, 2); 129359405d5eSPeter Brune *order = linesearch->order; 12943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129559405d5eSPeter Brune } 129659405d5eSPeter Brune 129759405d5eSPeter Brune /*@ 12981f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 129959405d5eSPeter Brune 130059405d5eSPeter Brune Input Parameters: 1301420bcc1bSBarry Smith + linesearch - the line search context 1302ceaaa498SBarry Smith - order - The order 130359405d5eSPeter Brune 130459405d5eSPeter Brune Level: intermediate 130559405d5eSPeter Brune 1306ceaaa498SBarry Smith Values for `order`\: 1307f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1308f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1309f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 131078bcb3b5SPeter Brune 1311420bcc1bSBarry Smith Options Database Key: 1312420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3) 1313420bcc1bSBarry Smith 1314ceaaa498SBarry Smith Note: 1315ceaaa498SBarry Smith These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP` 131659405d5eSPeter Brune 1317420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 131859405d5eSPeter Brune @*/ 1319d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1320d71ae5a4SJacob Faibussowitsch { 132159405d5eSPeter Brune PetscFunctionBegin; 132259405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 132359405d5eSPeter Brune linesearch->order = order; 13243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132559405d5eSPeter Brune } 132659405d5eSPeter Brune 1327f40b411bSPeter Brune /*@ 1328420bcc1bSBarry Smith SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1329f40b411bSPeter Brune 1330f6dfbefdSBarry Smith Not Collective 1331f6dfbefdSBarry Smith 1332f899ff85SJose E. Roman Input Parameter: 1333420bcc1bSBarry Smith . linesearch - the line search context 1334f40b411bSPeter Brune 1335f40b411bSPeter Brune Output Parameters: 1336f40b411bSPeter Brune + xnorm - The norm of the current solution 1337a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1338f40b411bSPeter Brune - ynorm - The norm of the current update 1339f40b411bSPeter Brune 134078bcb3b5SPeter Brune Level: developer 1341f40b411bSPeter Brune 1342420bcc1bSBarry Smith Notes: 1343420bcc1bSBarry Smith Some values may not be up-to-date at particular points in the code. 1344a68bbae5SBarry Smith 1345a68bbae5SBarry Smith This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1346a68bbae5SBarry Smith computed values. 1347a68bbae5SBarry Smith 1348420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1349f40b411bSPeter Brune @*/ 1350d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1351d71ae5a4SJacob Faibussowitsch { 1352bf7f4e0aSPeter Brune PetscFunctionBegin; 1353f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1354f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1355f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1356f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 13573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1358bf7f4e0aSPeter Brune } 1359bf7f4e0aSPeter Brune 1360f40b411bSPeter Brune /*@ 1361420bcc1bSBarry Smith SNESLineSearchSetNorms - Gets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1362f40b411bSPeter Brune 1363c3339decSBarry Smith Collective 1364f6dfbefdSBarry Smith 1365f40b411bSPeter Brune Input Parameters: 1366420bcc1bSBarry Smith + linesearch - the line search context 1367f40b411bSPeter Brune . xnorm - The norm of the current solution 1368a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1369f40b411bSPeter Brune - ynorm - The norm of the current update 1370f40b411bSPeter Brune 1371f6dfbefdSBarry Smith Level: developer 1372f40b411bSPeter Brune 1373420bcc1bSBarry Smith Note: 1374420bcc1bSBarry Smith This is called by the line search routines to store the values they have just computed 1375420bcc1bSBarry Smith 1376420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1377f40b411bSPeter Brune @*/ 1378d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1379d71ae5a4SJacob Faibussowitsch { 13806a388c36SPeter Brune PetscFunctionBegin; 1381f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 13826a388c36SPeter Brune linesearch->xnorm = xnorm; 13836a388c36SPeter Brune linesearch->fnorm = fnorm; 13846a388c36SPeter Brune linesearch->ynorm = ynorm; 13853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13866a388c36SPeter Brune } 13876a388c36SPeter Brune 1388f40b411bSPeter Brune /*@ 1389420bcc1bSBarry Smith SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`. 1390f40b411bSPeter Brune 13912fe279fdSBarry Smith Input Parameter: 1392420bcc1bSBarry Smith . linesearch - the line search context 1393f40b411bSPeter Brune 139420f4b53cSBarry Smith Options Database Key: 13958e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1396f40b411bSPeter Brune 1397f40b411bSPeter Brune Level: intermediate 1398f40b411bSPeter Brune 1399420bcc1bSBarry Smith Developer Note: 1400420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms 1401420bcc1bSBarry Smith 1402420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1403f40b411bSPeter Brune @*/ 1404d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1405d71ae5a4SJacob Faibussowitsch { 14069bd66eb0SPeter Brune SNES snes; 1407bf388a1fSBarry Smith 14086a388c36SPeter Brune PetscFunctionBegin; 14096a388c36SPeter Brune if (linesearch->norms) { 14109bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 14119566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 14129566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14139566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14149566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 14159bd66eb0SPeter Brune } else { 14169566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 14179566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14189566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14199566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 14209566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14219566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14226a388c36SPeter Brune } 14239bd66eb0SPeter Brune } 14243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14256a388c36SPeter Brune } 14266a388c36SPeter Brune 14276f263ca3SPeter Brune /*@ 14286f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 14296f263ca3SPeter Brune 14306f263ca3SPeter Brune Input Parameters: 1431420bcc1bSBarry Smith + linesearch - the line search context 143278bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 14336f263ca3SPeter Brune 143420f4b53cSBarry Smith Options Database Key: 1435420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search 14366f263ca3SPeter Brune 143720f4b53cSBarry Smith Level: intermediate 143820f4b53cSBarry Smith 1439f6dfbefdSBarry Smith Note: 1440f6dfbefdSBarry 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. 14416f263ca3SPeter Brune 1442420bcc1bSBarry Smith Developer Note: 1443420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms 1444420bcc1bSBarry Smith 1445420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 14466f263ca3SPeter Brune @*/ 1447d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1448d71ae5a4SJacob Faibussowitsch { 14496f263ca3SPeter Brune PetscFunctionBegin; 14506f263ca3SPeter Brune linesearch->norms = flg; 14513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14526f263ca3SPeter Brune } 14536f263ca3SPeter Brune 1454f40b411bSPeter Brune /*@ 1455f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1456f6dfbefdSBarry Smith 14578f14a041SBarry Smith Not Collective but the vectors are parallel 1458f40b411bSPeter Brune 1459f899ff85SJose E. Roman Input Parameter: 1460420bcc1bSBarry Smith . linesearch - the line search context 1461f40b411bSPeter Brune 1462f40b411bSPeter Brune Output Parameters: 14636232e825SPeter Brune + X - Solution vector 14646232e825SPeter Brune . F - Function vector 14656232e825SPeter Brune . Y - Search direction vector 14666232e825SPeter Brune . W - Solution work vector 14676232e825SPeter Brune - G - Function work vector 14686232e825SPeter Brune 146920f4b53cSBarry Smith Level: advanced 147020f4b53cSBarry Smith 14717bba9028SPeter Brune Notes: 147220f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a 147320f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the 147420f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the 14756232e825SPeter Brune function evaluated at the new solution. 1476f40b411bSPeter Brune 1477f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 14782a7a6963SBarry Smith 1479420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1480f40b411bSPeter Brune @*/ 1481d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1482d71ae5a4SJacob Faibussowitsch { 14836a388c36SPeter Brune PetscFunctionBegin; 1484f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14856a388c36SPeter Brune if (X) { 14864f572ea9SToby Isaac PetscAssertPointer(X, 2); 14876a388c36SPeter Brune *X = linesearch->vec_sol; 14886a388c36SPeter Brune } 14896a388c36SPeter Brune if (F) { 14904f572ea9SToby Isaac PetscAssertPointer(F, 3); 14916a388c36SPeter Brune *F = linesearch->vec_func; 14926a388c36SPeter Brune } 14936a388c36SPeter Brune if (Y) { 14944f572ea9SToby Isaac PetscAssertPointer(Y, 4); 14956a388c36SPeter Brune *Y = linesearch->vec_update; 14966a388c36SPeter Brune } 14976a388c36SPeter Brune if (W) { 14984f572ea9SToby Isaac PetscAssertPointer(W, 5); 14996a388c36SPeter Brune *W = linesearch->vec_sol_new; 15006a388c36SPeter Brune } 15016a388c36SPeter Brune if (G) { 15024f572ea9SToby Isaac PetscAssertPointer(G, 6); 15036a388c36SPeter Brune *G = linesearch->vec_func_new; 15046a388c36SPeter Brune } 15053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15066a388c36SPeter Brune } 15076a388c36SPeter Brune 1508f40b411bSPeter Brune /*@ 1509f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1510f6dfbefdSBarry Smith 1511c3339decSBarry Smith Logically Collective 1512f40b411bSPeter Brune 1513f40b411bSPeter Brune Input Parameters: 1514420bcc1bSBarry Smith + linesearch - the line search context 15156232e825SPeter Brune . X - Solution vector 15166232e825SPeter Brune . F - Function vector 15176232e825SPeter Brune . Y - Search direction vector 15186232e825SPeter Brune . W - Solution work vector 15196232e825SPeter Brune - G - Function work vector 1520f40b411bSPeter Brune 1521420bcc1bSBarry Smith Level: developer 1522f40b411bSPeter Brune 1523420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1524f40b411bSPeter Brune @*/ 1525d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1526d71ae5a4SJacob Faibussowitsch { 15276a388c36SPeter Brune PetscFunctionBegin; 1528f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15296a388c36SPeter Brune if (X) { 15306a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 15316a388c36SPeter Brune linesearch->vec_sol = X; 15326a388c36SPeter Brune } 15336a388c36SPeter Brune if (F) { 15346a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 15356a388c36SPeter Brune linesearch->vec_func = F; 15366a388c36SPeter Brune } 15376a388c36SPeter Brune if (Y) { 15386a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 15396a388c36SPeter Brune linesearch->vec_update = Y; 15406a388c36SPeter Brune } 15416a388c36SPeter Brune if (W) { 15426a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 15436a388c36SPeter Brune linesearch->vec_sol_new = W; 15446a388c36SPeter Brune } 15456a388c36SPeter Brune if (G) { 15466a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 15476a388c36SPeter Brune linesearch->vec_func_new = G; 15486a388c36SPeter Brune } 15493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15506a388c36SPeter Brune } 15516a388c36SPeter Brune 1552e7058c64SPeter Brune /*@C 1553f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1554f6dfbefdSBarry Smith `SNESLineSearch` options in the database. 1555e7058c64SPeter Brune 1556c3339decSBarry Smith Logically Collective 1557e7058c64SPeter Brune 1558e7058c64SPeter Brune Input Parameters: 1559f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1560e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1561e7058c64SPeter Brune 156220f4b53cSBarry Smith Level: advanced 156320f4b53cSBarry Smith 1564f6dfbefdSBarry Smith Note: 1565e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1566e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1567e7058c64SPeter Brune 1568420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1569e7058c64SPeter Brune @*/ 1570d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1571d71ae5a4SJacob Faibussowitsch { 1572e7058c64SPeter Brune PetscFunctionBegin; 1573f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15749566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 15753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1576e7058c64SPeter Brune } 1577e7058c64SPeter Brune 1578e7058c64SPeter Brune /*@C 1579f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1580f1c6b773SPeter Brune SNESLineSearch options in the database. 1581e7058c64SPeter Brune 1582e7058c64SPeter Brune Not Collective 1583e7058c64SPeter Brune 1584e7058c64SPeter Brune Input Parameter: 1585f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 1586e7058c64SPeter Brune 1587e7058c64SPeter Brune Output Parameter: 1588e7058c64SPeter Brune . prefix - pointer to the prefix string used 1589e7058c64SPeter Brune 1590e7058c64SPeter Brune Level: advanced 1591e7058c64SPeter Brune 1592e4094ef1SJacob Faibussowitsch Fortran Notes: 159320f4b53cSBarry Smith The user should pass in a string 'prefix' of 159420f4b53cSBarry Smith sufficient length to hold the prefix. 159520f4b53cSBarry Smith 1596420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1597e7058c64SPeter Brune @*/ 1598d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1599d71ae5a4SJacob Faibussowitsch { 1600e7058c64SPeter Brune PetscFunctionBegin; 1601f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 16033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1604e7058c64SPeter Brune } 1605bf7f4e0aSPeter Brune 16068d359177SBarry Smith /*@C 1607f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1608f40b411bSPeter Brune 1609d8d19677SJose E. Roman Input Parameters: 1610f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1611f40b411bSPeter Brune - nwork - the number of work vectors 1612f40b411bSPeter Brune 1613f40b411bSPeter Brune Level: developer 1614f40b411bSPeter Brune 1615420bcc1bSBarry Smith Developer Note: 1616420bcc1bSBarry Smith This is called from within the set up routines for each of the line search types `SNESLineSearchType` 1617420bcc1bSBarry Smith 1618420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()` 1619f40b411bSPeter Brune @*/ 1620d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1621d71ae5a4SJacob Faibussowitsch { 1622bf7f4e0aSPeter Brune PetscFunctionBegin; 16230fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 16249566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 16253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1626bf7f4e0aSPeter Brune } 1627bf7f4e0aSPeter Brune 1628f40b411bSPeter Brune /*@ 1629422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1630f40b411bSPeter Brune 16312fe279fdSBarry Smith Input Parameter: 1632420bcc1bSBarry Smith . linesearch - the line search context 1633f40b411bSPeter Brune 16342fe279fdSBarry Smith Output Parameter: 1635422a814eSBarry Smith . result - The success or failure status 1636f40b411bSPeter Brune 163720f4b53cSBarry Smith Level: developer 163820f4b53cSBarry Smith 1639f6dfbefdSBarry Smith Note: 1640420bcc1bSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed 1641420bcc1bSBarry Smith (and set into the `SNES` convergence accordingly). 1642cd7522eaSPeter Brune 1643420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1644f40b411bSPeter Brune @*/ 1645d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1646d71ae5a4SJacob Faibussowitsch { 1647bf7f4e0aSPeter Brune PetscFunctionBegin; 1648f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16494f572ea9SToby Isaac PetscAssertPointer(result, 2); 1650422a814eSBarry Smith *result = linesearch->result; 16513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1652bf7f4e0aSPeter Brune } 1653bf7f4e0aSPeter Brune 1654f40b411bSPeter Brune /*@ 1655420bcc1bSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the line search application 1656f40b411bSPeter Brune 1657f40b411bSPeter Brune Input Parameters: 1658420bcc1bSBarry Smith + linesearch - the line search context 1659422a814eSBarry Smith - result - The success or failure status 1660f40b411bSPeter Brune 166120f4b53cSBarry Smith Level: developer 166220f4b53cSBarry Smith 1663f6dfbefdSBarry Smith Note: 1664420bcc1bSBarry Smith This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1665cd7522eaSPeter Brune the success or failure of the line search method. 1666cd7522eaSPeter Brune 1667420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1668f40b411bSPeter Brune @*/ 1669d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1670d71ae5a4SJacob Faibussowitsch { 16716a388c36SPeter Brune PetscFunctionBegin; 16725fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1673422a814eSBarry Smith linesearch->result = result; 16743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16756a388c36SPeter Brune } 16766a388c36SPeter Brune 1677ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 16789bd66eb0SPeter Brune /*@C 1679f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16809bd66eb0SPeter Brune 1681c3339decSBarry Smith Logically Collective 1682f6dfbefdSBarry Smith 16839bd66eb0SPeter Brune Input Parameters: 1684ceaaa498SBarry Smith + linesearch - the linesearch object 1685420bcc1bSBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchSetVIFunctions` for calling sequence 1686420bcc1bSBarry Smith - normfunc - function for computing the norm of an active set, see `SNESLineSearchSetVIFunctions ` for calling sequence 16879bd66eb0SPeter Brune 1688f6dfbefdSBarry Smith Level: advanced 16899bd66eb0SPeter Brune 1690ceaaa498SBarry Smith Notes: 169120f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 169220f4b53cSBarry Smith 169320f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated 169420f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`. 169520f4b53cSBarry Smith 1696420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 16979bd66eb0SPeter Brune @*/ 1698d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1699d71ae5a4SJacob Faibussowitsch { 17009bd66eb0SPeter Brune PetscFunctionBegin; 1701f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 17029bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 17039bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 17043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17059bd66eb0SPeter Brune } 17069bd66eb0SPeter Brune 17079bd66eb0SPeter Brune /*@C 1708f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 17099bd66eb0SPeter Brune 1710f6dfbefdSBarry Smith Not Collective 1711f6dfbefdSBarry Smith 1712f899ff85SJose E. Roman Input Parameter: 1713f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()` 17149bd66eb0SPeter Brune 17159bd66eb0SPeter Brune Output Parameters: 17169bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 17179bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 17189bd66eb0SPeter Brune 1719f6dfbefdSBarry Smith Level: advanced 17209bd66eb0SPeter Brune 1721420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 17229bd66eb0SPeter Brune @*/ 1723d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1724d71ae5a4SJacob Faibussowitsch { 17259bd66eb0SPeter Brune PetscFunctionBegin; 17269bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 17279bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 17283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17299bd66eb0SPeter Brune } 17309bd66eb0SPeter Brune 1731bf7f4e0aSPeter Brune /*@C 1732420bcc1bSBarry Smith SNESLineSearchRegister - register a line search type `SNESLineSearchType` 1733ceaaa498SBarry Smith 1734ceaaa498SBarry Smith Input Parameters: 1735420bcc1bSBarry Smith + sname - name of the `SNESLineSearchType()` 1736ceaaa498SBarry Smith - function - the creation function for that type 1737ceaaa498SBarry Smith 1738ceaaa498SBarry Smith Calling sequence of `function`: 1739ceaaa498SBarry Smith . ls - the line search context 1740bf7f4e0aSPeter Brune 1741bf7f4e0aSPeter Brune Level: advanced 1742f6dfbefdSBarry Smith 1743420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1744bf7f4e0aSPeter Brune @*/ 1745ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls)) 1746d71ae5a4SJacob Faibussowitsch { 1747bf7f4e0aSPeter Brune PetscFunctionBegin; 17489566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 17499566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 17503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1751bf7f4e0aSPeter Brune } 1752