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) 8149abdd8aSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence 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 8720f4b53cSBarry Smith Level: intermediate 88dcf2fd19SBarry Smith 89f6dfbefdSBarry Smith Note: 90dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 91f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` multiple times; all will be called in the 92dcf2fd19SBarry Smith order in which they were set. 93dcf2fd19SBarry Smith 94420bcc1bSBarry Smith Fortran Note: 95f6dfbefdSBarry Smith Only a single monitor function can be set for each `SNESLineSearch` object 96dcf2fd19SBarry Smith 9749abdd8aSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`, `PetscCtxDestroyFn` 98dcf2fd19SBarry Smith @*/ 9949abdd8aSBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscCtxDestroyFn *monitordestroy) 100d71ae5a4SJacob Faibussowitsch { 10178064530SBarry Smith PetscInt i; 10278064530SBarry Smith PetscBool identical; 10378064530SBarry Smith 104dcf2fd19SBarry Smith PetscFunctionBegin; 105dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 10678064530SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 1079566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f, mctx, monitordestroy, (PetscErrorCode (*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical)); 1083ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 10978064530SBarry Smith } 1105f80ce2aSJacob Faibussowitsch PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 111dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 112dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 113835f2295SStefano Zampini ls->monitorcontext[ls->numbermonitors++] = mctx; 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115dcf2fd19SBarry Smith } 116dcf2fd19SBarry Smith 117dcf2fd19SBarry Smith /*@C 118f6dfbefdSBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries 119dcf2fd19SBarry Smith 120c3339decSBarry Smith Collective 121dcf2fd19SBarry Smith 122dcf2fd19SBarry Smith Input Parameters: 123420bcc1bSBarry Smith + ls - the `SNESLineSearch` object 124f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat` 125dcf2fd19SBarry Smith 126420bcc1bSBarry Smith Options Database Key: 127420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 128dcf2fd19SBarry Smith 129420bcc1bSBarry Smith Level: developer 130420bcc1bSBarry Smith 131420bcc1bSBarry Smith This is not normally called directly but is passed to `SNESLineSearchMonitorSet()` 132420bcc1bSBarry Smith 133420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()` 134dcf2fd19SBarry Smith @*/ 135d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) 136d71ae5a4SJacob Faibussowitsch { 137d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 138dcf2fd19SBarry Smith Vec Y, W, G; 139dcf2fd19SBarry Smith 140dcf2fd19SBarry Smith PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n")); 1449566063dSJacob Faibussowitsch PetscCall(VecView(Y, viewer)); 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n")); 1469566063dSJacob Faibussowitsch PetscCall(VecView(W, viewer)); 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n")); 1489566063dSJacob Faibussowitsch PetscCall(VecView(G, viewer)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151dcf2fd19SBarry Smith } 152dcf2fd19SBarry Smith 153f40b411bSPeter Brune /*@ 154420bcc1bSBarry Smith SNESLineSearchCreate - Creates a `SNESLineSearch` context. 155f40b411bSPeter Brune 156f6dfbefdSBarry Smith Logically Collective 157f40b411bSPeter Brune 1582fe279fdSBarry Smith Input Parameter: 159f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context). 160f40b411bSPeter Brune 1612fe279fdSBarry Smith Output Parameter: 1628e557f58SPeter Brune . outlinesearch - the new line search context 163f40b411bSPeter Brune 164162e0bf5SPeter Brune Level: developer 165162e0bf5SPeter Brune 166f6dfbefdSBarry Smith Note: 167420bcc1bSBarry Smith The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance 168f6dfbefdSBarry Smith already associated with the `SNES`. 169f40b411bSPeter Brune 170420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()` 171f40b411bSPeter Brune @*/ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 173d71ae5a4SJacob Faibussowitsch { 174f1c6b773SPeter Brune SNESLineSearch linesearch; 175bf388a1fSBarry Smith 176bf7f4e0aSPeter Brune PetscFunctionBegin; 1774f572ea9SToby Isaac PetscAssertPointer(outlinesearch, 2); 1789566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 179f5af7f23SKarl Rupp 1809566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView)); 1810298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1820298fd71SBarry Smith linesearch->vec_func_new = NULL; 1830298fd71SBarry Smith linesearch->vec_sol = NULL; 1840298fd71SBarry Smith linesearch->vec_func = NULL; 1850298fd71SBarry Smith linesearch->vec_update = NULL; 1869bd66eb0SPeter Brune 187bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 188bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 189bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 190bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 191422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 192bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 193bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 194bf7f4e0aSPeter Brune linesearch->damping = 1.0; 195a99ef635SJonas Heinzmann linesearch->maxlambda = 1.0; 196a99ef635SJonas Heinzmann linesearch->minlambda = 1e-12; 197516fe3c3SPeter Brune linesearch->rtol = 1e-8; 198516fe3c3SPeter Brune linesearch->atol = 1e-15; 199516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2000298fd71SBarry Smith linesearch->precheckctx = NULL; 2010298fd71SBarry Smith linesearch->postcheckctx = NULL; 202a99ef635SJonas Heinzmann linesearch->max_it = 1; 203bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2043add74b1SBarry Smith linesearch->monitor = NULL; 205bf7f4e0aSPeter Brune *outlinesearch = linesearch; 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207bf7f4e0aSPeter Brune } 208bf7f4e0aSPeter Brune 209f40b411bSPeter Brune /*@ 21078bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 21178bcb3b5SPeter Brune any required vectors. 212f40b411bSPeter Brune 213c3339decSBarry Smith Collective 214f40b411bSPeter Brune 2152fe279fdSBarry Smith Input Parameter: 216f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 217f40b411bSPeter Brune 21820f4b53cSBarry Smith Level: advanced 21920f4b53cSBarry Smith 220f6dfbefdSBarry Smith Note: 221f6dfbefdSBarry Smith For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`. 222cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 22378bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 224f6dfbefdSBarry Smith of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be 225cd7522eaSPeter Brune allocated upfront. 226cd7522eaSPeter Brune 227420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()` 228f40b411bSPeter Brune @*/ 229d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 230d71ae5a4SJacob Faibussowitsch { 231bf7f4e0aSPeter Brune PetscFunctionBegin; 23248a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 233bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 23448a46eb9SPierre Jolivet if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 23548a46eb9SPierre Jolivet if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 236213acdd3SPierre Jolivet PetscTryTypeMethod(linesearch, setup); 2379566063dSJacob Faibussowitsch if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction)); 238bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 239bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 240bf7f4e0aSPeter Brune } 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242bf7f4e0aSPeter Brune } 243bf7f4e0aSPeter Brune 244f40b411bSPeter Brune /*@ 245f6dfbefdSBarry Smith SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search. 246f40b411bSPeter Brune 247c3339decSBarry Smith Collective 248f40b411bSPeter Brune 2492fe279fdSBarry Smith Input Parameter: 250f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 251f40b411bSPeter Brune 25220f4b53cSBarry Smith Level: developer 25320f4b53cSBarry Smith 254f6dfbefdSBarry Smith Note: 255f6dfbefdSBarry Smith Usually only called by `SNESReset()` 256f190f2fcSBarry Smith 257420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 258f40b411bSPeter Brune @*/ 259d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 260d71ae5a4SJacob Faibussowitsch { 261bf7f4e0aSPeter Brune PetscFunctionBegin; 262213acdd3SPierre Jolivet PetscTryTypeMethod(linesearch, reset); 263f5af7f23SKarl Rupp 2649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_sol_new)); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_func_new)); 266bf7f4e0aSPeter Brune 2679566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 268f5af7f23SKarl Rupp 269bf7f4e0aSPeter Brune linesearch->nwork = 0; 270bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 272bf7f4e0aSPeter Brune } 273bf7f4e0aSPeter Brune 274ed07d7d7SPeter Brune /*@C 275f6dfbefdSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search 276f6dfbefdSBarry Smith ` 277e4094ef1SJacob Faibussowitsch 278ed07d7d7SPeter Brune Input Parameters: 279e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context 280e4094ef1SJacob Faibussowitsch - func - function evaluation routine, this is usually the function provided with `SNESSetFunction()` 281ed07d7d7SPeter Brune 282420bcc1bSBarry Smith Calling sequence of `func`: 283420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with 284420bcc1bSBarry Smith . x - the input vector 285420bcc1bSBarry Smith - f - the computed value of the function 286420bcc1bSBarry Smith 287ed07d7d7SPeter Brune Level: developer 288ed07d7d7SPeter Brune 289420bcc1bSBarry Smith Note: 290420bcc1bSBarry Smith By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed 291420bcc1bSBarry Smith 292420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()` 293ed07d7d7SPeter Brune @*/ 294420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f)) 295d71ae5a4SJacob Faibussowitsch { 296ed07d7d7SPeter Brune PetscFunctionBegin; 297ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 298ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 2993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 300ed07d7d7SPeter Brune } 301ed07d7d7SPeter Brune 30286d74e61SPeter Brune /*@C 303420bcc1bSBarry Smith SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but 304420bcc1bSBarry Smith before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that 305f190f2fcSBarry Smith determined the search direction. 30686d74e61SPeter Brune 307c3339decSBarry Smith Logically Collective 30886d74e61SPeter Brune 30986d74e61SPeter Brune Input Parameters: 310f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 311420bcc1bSBarry Smith . func - [optional] function evaluation routine 31220f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 31386d74e61SPeter Brune 314420bcc1bSBarry Smith Calling sequence of `func`: 315420bcc1bSBarry Smith + ls - the `SNESLineSearch` context 316420bcc1bSBarry Smith . x - the current solution 317a678f235SPierre Jolivet . d - the current search direction 318420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed 319420bcc1bSBarry Smith - ctx - the context passed to `SNESLineSearchSetPreCheck()` 320420bcc1bSBarry Smith 32186d74e61SPeter Brune Level: intermediate 32286d74e61SPeter Brune 323f6dfbefdSBarry Smith Note: 324420bcc1bSBarry Smith Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete. 325f0b84518SBarry Smith 326f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 327f0b84518SBarry Smith 328420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 329869ba2dcSBarry Smith `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 330f0b84518SBarry Smith 33186d74e61SPeter Brune @*/ 332420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx) 333d71ae5a4SJacob Faibussowitsch { 3349bd66eb0SPeter Brune PetscFunctionBegin; 335f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 336f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 33786d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33986d74e61SPeter Brune } 34086d74e61SPeter Brune 34186d74e61SPeter Brune /*@C 34253deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 34386d74e61SPeter Brune 344f899ff85SJose E. Roman Input Parameter: 345f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 34686d74e61SPeter Brune 34786d74e61SPeter Brune Output Parameters: 348420bcc1bSBarry Smith + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchSetPreCheck()` 34920f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 35086d74e61SPeter Brune 35186d74e61SPeter Brune Level: intermediate 35286d74e61SPeter Brune 353420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 35486d74e61SPeter Brune @*/ 355d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) 356d71ae5a4SJacob Faibussowitsch { 3579bd66eb0SPeter Brune PetscFunctionBegin; 358f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 359f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 36086d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36286d74e61SPeter Brune } 36386d74e61SPeter Brune 36486d74e61SPeter Brune /*@C 365f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 366f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 36786d74e61SPeter Brune 36820f4b53cSBarry Smith Logically Collective 36986d74e61SPeter Brune 37086d74e61SPeter Brune Input Parameters: 371f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 372420bcc1bSBarry Smith . func - [optional] function evaluation routine 37320f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 37486d74e61SPeter Brune 375420bcc1bSBarry Smith Calling sequence of `func`: 376420bcc1bSBarry Smith + ls - the `SNESLineSearch` context 377420bcc1bSBarry Smith . x - the current solution 378a678f235SPierre Jolivet . d - the current search direction 379a678f235SPierre Jolivet . w - $ w = x + lambda*d $ for some lambda 380420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed 381420bcc1bSBarry Smith . changed_w - indicates `w` has been changed 382420bcc1bSBarry Smith - ctx - the context passed to `SNESLineSearchSetPreCheck()` 383420bcc1bSBarry Smith 38486d74e61SPeter Brune Level: intermediate 38586d74e61SPeter Brune 386f0b84518SBarry Smith Notes: 3876b095a85SStefano Zampini Use `SNESLineSearchSetPreCheck()` to change the step before the line search is completed. 3886b095a85SStefano Zampini The calling sequence of the callback does not contain the current scaling factor. To access the value, use `SNESLineSearchGetLambda()`. 389f0b84518SBarry Smith 390f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 391f0b84518SBarry Smith 392420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`, 393e4094ef1SJacob Faibussowitsch `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 39486d74e61SPeter Brune @*/ 395420bcc1bSBarry 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) 396d71ae5a4SJacob Faibussowitsch { 39786d74e61SPeter Brune PetscFunctionBegin; 398f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 399f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 40086d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40286d74e61SPeter Brune } 40386d74e61SPeter Brune 40486d74e61SPeter Brune /*@C 405f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 40686d74e61SPeter Brune 407f899ff85SJose E. Roman Input Parameter: 408f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 40986d74e61SPeter Brune 41086d74e61SPeter Brune Output Parameters: 411420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()` 41220f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 41386d74e61SPeter Brune 41486d74e61SPeter Brune Level: intermediate 41586d74e61SPeter Brune 416420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 41786d74e61SPeter Brune @*/ 418d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 419d71ae5a4SJacob Faibussowitsch { 4209bd66eb0SPeter Brune PetscFunctionBegin; 421f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 422f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 42386d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 4243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42586d74e61SPeter Brune } 42686d74e61SPeter Brune 427f40b411bSPeter Brune /*@ 428f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 429f40b411bSPeter Brune 430c3339decSBarry Smith Logically Collective 431f40b411bSPeter Brune 432f40b411bSPeter Brune Input Parameters: 4337b1df9c1SPeter Brune + linesearch - The linesearch instance. 4347b1df9c1SPeter Brune . X - The current solution 4357b1df9c1SPeter Brune - Y - The step direction 436f40b411bSPeter Brune 4372fe279fdSBarry Smith Output Parameter: 438420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y` 439f40b411bSPeter Brune 440f0b84518SBarry Smith Level: advanced 441f40b411bSPeter Brune 442f6dfbefdSBarry Smith Note: 443420bcc1bSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines 444f6dfbefdSBarry Smith 445420bcc1bSBarry Smith Developer Note: 446420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided 447420bcc1bSBarry Smith 448420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, 449fe8e7dddSPierre Jolivet `SNESLineSearchGetPostCheck()` 450f40b411bSPeter Brune @*/ 451d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) 452d71ae5a4SJacob Faibussowitsch { 453bf7f4e0aSPeter Brune PetscFunctionBegin; 454bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 4556b2b7091SBarry Smith if (linesearch->ops->precheck) { 456dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx); 45738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed, 4); 458bf7f4e0aSPeter Brune } 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 460bf7f4e0aSPeter Brune } 461bf7f4e0aSPeter Brune 462f40b411bSPeter Brune /*@ 463ef46b1a6SStefano Zampini SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 464f40b411bSPeter Brune 465c3339decSBarry Smith Logically Collective 466f40b411bSPeter Brune 467f40b411bSPeter Brune Input Parameters: 4687b1df9c1SPeter Brune + linesearch - The line search context 4697b1df9c1SPeter Brune . X - The last solution 4707b1df9c1SPeter Brune . Y - The step direction 4716b095a85SStefano Zampini - W - The updated solution, `W = X - lambda * Y` for some lambda 472f40b411bSPeter Brune 473f40b411bSPeter Brune Output Parameters: 4746b095a85SStefano Zampini + changed_Y - Indicator if the direction `Y` has been changed. 475420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed. 476f40b411bSPeter Brune 47720f4b53cSBarry Smith Level: developer 47820f4b53cSBarry Smith 479f6dfbefdSBarry Smith Note: 480420bcc1bSBarry Smith This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines 481f6dfbefdSBarry Smith 482420bcc1bSBarry Smith Developer Note: 483420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided 484420bcc1bSBarry Smith 485420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 486f40b411bSPeter Brune @*/ 487d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 488d71ae5a4SJacob Faibussowitsch { 489bf7f4e0aSPeter Brune PetscFunctionBegin; 490bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 491bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4926b2b7091SBarry Smith if (linesearch->ops->postcheck) { 493dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx); 49438bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5); 49538bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6); 49686d74e61SPeter Brune } 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49886d74e61SPeter Brune } 49986d74e61SPeter Brune 50086d74e61SPeter Brune /*@C 5011d27aa22SBarry Smith SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time` 50286d74e61SPeter Brune 503c3339decSBarry Smith Logically Collective 50486d74e61SPeter Brune 5054165533cSJose E. Roman Input Parameters: 506420bcc1bSBarry Smith + linesearch - the line search context 50786d74e61SPeter Brune . X - base state for this step 508907376e6SBarry Smith - ctx - context for this function 50986d74e61SPeter Brune 51097bb3fdcSJose E. Roman Input/Output Parameter: 51197bb3fdcSJose E. Roman . Y - correction, possibly modified 51297bb3fdcSJose E. Roman 51397bb3fdcSJose E. Roman Output Parameter: 514420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified 51586d74e61SPeter Brune 516420bcc1bSBarry Smith Options Database Keys: 517cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 518cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 51986d74e61SPeter Brune 52086d74e61SPeter Brune Level: advanced 52186d74e61SPeter Brune 52286d74e61SPeter Brune Notes: 523f6dfbefdSBarry Smith This function should be passed to `SNESLineSearchSetPreCheck()` 52486d74e61SPeter Brune 52586d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 5261d27aa22SBarry Smith so the Picard linearization should be provided in place of the "Jacobian" {cite}`hindmarsh1996time`. This correction 52786d74e61SPeter Brune is generally not useful when using a Newton linearization. 52886d74e61SPeter Brune 529420bcc1bSBarry Smith Developer Note: 530420bcc1bSBarry Smith The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided 531420bcc1bSBarry Smith 532420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 53386d74e61SPeter Brune @*/ 534d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) 535d71ae5a4SJacob Faibussowitsch { 53686d74e61SPeter Brune PetscReal angle = *(PetscReal *)linesearch->precheckctx; 53786d74e61SPeter Brune Vec Ylast; 53886d74e61SPeter Brune PetscScalar dot; 53986d74e61SPeter Brune PetscInt iter; 54086d74e61SPeter Brune PetscReal ynorm, ylastnorm, theta, angle_radians; 54186d74e61SPeter Brune SNES snes; 54286d74e61SPeter Brune 54386d74e61SPeter Brune PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast)); 54686d74e61SPeter Brune if (!Ylast) { 5479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Ylast)); 5489566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast)); 5499566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Ylast)); 55086d74e61SPeter Brune } 5519566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 55286d74e61SPeter Brune if (iter < 2) { 5539566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 55486d74e61SPeter Brune *changed = PETSC_FALSE; 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55686d74e61SPeter Brune } 55786d74e61SPeter Brune 5589566063dSJacob Faibussowitsch PetscCall(VecDot(Y, Ylast, &dot)); 5599566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); 5609566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 561c197fec0SMatthew Knepley if (ynorm == 0. || ylastnorm == 0.) { 562c197fec0SMatthew Knepley *changed = PETSC_FALSE; 563c197fec0SMatthew Knepley PetscFunctionReturn(PETSC_SUCCESS); 564c197fec0SMatthew Knepley } 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)); 576835f2295SStefano Zampini 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 { 579835f2295SStefano Zampini 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 5916b095a85SStefano Zampini Input Parameter: 5926b095a85SStefano Zampini . linesearch - The line search context 593f40b411bSPeter Brune 5946b867d5aSJose E. Roman Input/Output Parameters: 5956b867d5aSJose E. Roman + X - The current solution, on output the new solution 596420bcc1bSBarry Smith . F - The current function value, on output the new function value at the solution value `X` 5971bd63e3eSSatish Balay . fnorm - The current norm of `F`, on output the new norm of `F` 598a99ef635SJonas Heinzmann - Y - The current search direction, on output the direction determined by the linesearch, i.e. `Xnew = Xold - lambda*Y` 599f40b411bSPeter Brune 600cd7522eaSPeter Brune Options Database Keys: 601a99ef635SJonas Heinzmann + -snes_linesearch_type - basic (or equivalently none), bt, secant, cp, nleqerr, bisection, 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()) 605a99ef635SJonas Heinzmann . -snes_linesearch_keeplambda - Keep the previous `lambda` 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 614a99ef635SJonas Heinzmann an optimal damping parameter (that is `lambda`) 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 6181bd63e3eSSatish Balay .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `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); 671f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1); 672f4f49eeaSPierre Jolivet 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); 679648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 680f4f49eeaSPierre Jolivet 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; 710648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&linesearch->monitor)); 711dcf2fd19SBarry Smith linesearch->monitor = viewer; 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 713bf7f4e0aSPeter Brune } 714bf7f4e0aSPeter Brune 715f40b411bSPeter Brune /*@ 716420bcc1bSBarry Smith SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()` 717f6dfbefdSBarry Smith 718c3339decSBarry Smith Logically Collective 7196a388c36SPeter Brune 720f190f2fcSBarry Smith Input Parameter: 721420bcc1bSBarry Smith . linesearch - the line search context 722f40b411bSPeter Brune 723f190f2fcSBarry Smith Output Parameter: 7248e557f58SPeter Brune . monitor - monitor context 725f40b411bSPeter Brune 726f40b411bSPeter Brune Level: intermediate 727f40b411bSPeter Brune 728420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 729f40b411bSPeter Brune @*/ 730d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 731d71ae5a4SJacob Faibussowitsch { 7326a388c36SPeter Brune PetscFunctionBegin; 733f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 7346a388c36SPeter Brune *monitor = linesearch->monitor; 7353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7366a388c36SPeter Brune } 7376a388c36SPeter Brune 738dcf2fd19SBarry Smith /*@C 739420bcc1bSBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database 740dcf2fd19SBarry Smith 741c3339decSBarry Smith Collective 742dcf2fd19SBarry Smith 743dcf2fd19SBarry Smith Input Parameters: 744420bcc1bSBarry Smith + ls - `SNESLineSearch` object to monitor 745420bcc1bSBarry Smith . name - the monitor type 746dcf2fd19SBarry Smith . help - message indicating what monitoring is done 747dcf2fd19SBarry Smith . manual - manual page for the monitor 74849abdd8aSBarry Smith . monitor - the monitor function, must use `PetscViewerAndFormat` as its context 749f6dfbefdSBarry 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` 750dcf2fd19SBarry Smith 751420bcc1bSBarry Smith Calling sequence of `monitor`: 752420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored 753420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used 754dcf2fd19SBarry Smith 755420bcc1bSBarry Smith Calling sequence of `monitorsetup`: 756420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored 757420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used 758420bcc1bSBarry Smith 759420bcc1bSBarry Smith Level: advanced 760420bcc1bSBarry Smith 761648c30bcSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 762db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 763e4094ef1SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 764db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 765c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 766db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 767db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 768dcf2fd19SBarry Smith @*/ 769420bcc1bSBarry 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)) 770d71ae5a4SJacob Faibussowitsch { 771dcf2fd19SBarry Smith PetscViewer viewer; 772dcf2fd19SBarry Smith PetscViewerFormat format; 773dcf2fd19SBarry Smith PetscBool flg; 774dcf2fd19SBarry Smith 775dcf2fd19SBarry Smith PetscFunctionBegin; 776648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg)); 777dcf2fd19SBarry Smith if (flg) { 778d12e167eSBarry Smith PetscViewerAndFormat *vf; 7799566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 780648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 7811baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls, vf)); 78249abdd8aSBarry Smith PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode (*)(SNESLineSearch, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy)); 783dcf2fd19SBarry Smith } 7843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 785dcf2fd19SBarry Smith } 786dcf2fd19SBarry Smith 787f40b411bSPeter Brune /*@ 788f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 789f40b411bSPeter Brune 790c3339decSBarry Smith Logically Collective 791f6dfbefdSBarry Smith 792f899ff85SJose E. Roman Input Parameter: 793420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context 794f40b411bSPeter Brune 795f40b411bSPeter Brune Options Database Keys: 796a99ef635SJonas Heinzmann + -snes_linesearch_type <type> - basic (or equivalently none), `bt`, `secant`, `cp`, `nleqerr`, `bisection`, `shell` 797a99ef635SJonas Heinzmann . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`bt` supports 1, 2 or 3) 798f6dfbefdSBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`) 799a99ef635SJonas Heinzmann . -snes_linesearch_minlambda - The minimum `lambda` 800a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda - The maximum `lambda` 8011a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 8021a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 803a99ef635SJonas Heinzmann . -snes_linesearch_ltol - Change in `lambda` tolerance for iterative line searches 8041a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 805dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 806dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 8078e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 808a99ef635SJonas Heinzmann . -snes_linesearch_keeplambda - Keep the previous `lambda` as the initial guess. 8091a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 810d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 811f40b411bSPeter Brune 812f40b411bSPeter Brune Level: intermediate 813f40b411bSPeter Brune 81462842358SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 815db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 816f40b411bSPeter Brune @*/ 817d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 818d71ae5a4SJacob Faibussowitsch { 8191a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 820bf7f4e0aSPeter Brune char type[256]; 821bf7f4e0aSPeter Brune PetscBool flg, set; 822dcf2fd19SBarry Smith PetscViewer viewer; 823bf388a1fSBarry Smith 824bf7f4e0aSPeter Brune PetscFunctionBegin; 8259566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll()); 826bf7f4e0aSPeter Brune 827d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch); 828f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 8299566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg)); 830bf7f4e0aSPeter Brune if (flg) { 8319566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, type)); 832bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 8339566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, deft)); 834bf7f4e0aSPeter Brune } 835bf7f4e0aSPeter Brune 836648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set)); 837648c30bcSBarry Smith if (set) PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer)); 8389566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL)); 839bf7f4e0aSPeter Brune 8401a9b3a06SPeter Brune /* tolerances */ 841a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum lambda", "SNESLineSearchSetTolerances", linesearch->minlambda, &linesearch->minlambda, NULL)); 842a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_maxlambda", "Maximum lambda", "SNESLineSearchSetTolerances", linesearch->maxlambda, &linesearch->maxlambda, NULL)); 8439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL)); 8449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL)); 8459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL)); 846a99ef635SJonas Heinzmann PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_it, &linesearch->max_it, NULL)); 847a99ef635SJonas Heinzmann 848a99ef635SJonas Heinzmann /* deprecated options */ 849a99ef635SJonas Heinzmann PetscCall(PetscOptionsDeprecated("-snes_linesearch_maxstep", "-snes_linesearch_maxlambda", "3.24.0", NULL)); 8507a35526eSPeter Brune 8511a9b3a06SPeter Brune /* damping parameters */ 852a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping (and depending on chosen line search initial lambda guess)", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 8531a9b3a06SPeter Brune 8549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 8551a9b3a06SPeter Brune 8561a9b3a06SPeter Brune /* precheck */ 8579566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 8587a35526eSPeter Brune if (set) { 8597a35526eSPeter Brune if (flg) { 8607a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 861f5af7f23SKarl Rupp 862d0609cedSBarry 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)); 8639566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 8647a35526eSPeter Brune } else { 8659566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 8667a35526eSPeter Brune } 8677a35526eSPeter Brune } 8689566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 8699566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 8707a35526eSPeter Brune 871dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 8725a9b6599SPeter Brune 873dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 874d0609cedSBarry Smith PetscOptionsEnd(); 8753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 876bf7f4e0aSPeter Brune } 877bf7f4e0aSPeter Brune 878f40b411bSPeter Brune /*@ 879f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 880f40b411bSPeter Brune 88120f4b53cSBarry Smith Logically Collective 88220f4b53cSBarry Smith 883f40b411bSPeter Brune Input Parameters: 8842fe279fdSBarry Smith + linesearch - line search context 885420bcc1bSBarry Smith - viewer - the `PetscViewer` to display the line search information to 886f40b411bSPeter Brune 887f40b411bSPeter Brune Level: intermediate 888f40b411bSPeter Brune 889420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 890f40b411bSPeter Brune @*/ 891d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 892d71ae5a4SJacob Faibussowitsch { 893*9f196a02SMartin Diehl PetscBool isascii; 894bf388a1fSBarry Smith 895bf7f4e0aSPeter Brune PetscFunctionBegin; 8967f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 89748a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 8987f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 8997f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2); 900f40b411bSPeter Brune 901*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 902*9f196a02SMartin Diehl if (isascii) { 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 9049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 905dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer); 9069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 907a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(viewer, " maxlambda=%e, minlambda=%e\n", (double)linesearch->maxlambda, (double)linesearch->minlambda)); 9089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 909a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_it)); 9106b2b7091SBarry Smith if (linesearch->ops->precheck) { 9116b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 91263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 9137f1410a3SPeter Brune } else { 91463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 9157f1410a3SPeter Brune } 9167f1410a3SPeter Brune } 91748a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 9187f1410a3SPeter Brune } 9193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 920bf7f4e0aSPeter Brune } 921bf7f4e0aSPeter Brune 922cc4c1da9SBarry Smith /*@ 923420bcc1bSBarry Smith SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch` 924a80ff896SJed Brown 925c3339decSBarry Smith Logically Collective 926a80ff896SJed Brown 9272fe279fdSBarry Smith Input Parameter: 928420bcc1bSBarry Smith . linesearch - the line search context 929a80ff896SJed Brown 9302fe279fdSBarry Smith Output Parameter: 9312fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set 932a80ff896SJed Brown 933a80ff896SJed Brown Level: intermediate 934a80ff896SJed Brown 935420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 936a80ff896SJed Brown @*/ 937d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 938d71ae5a4SJacob Faibussowitsch { 939a80ff896SJed Brown PetscFunctionBegin; 940a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9414f572ea9SToby Isaac PetscAssertPointer(type, 2); 942a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 9433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 944a80ff896SJed Brown } 945a80ff896SJed Brown 946cc4c1da9SBarry Smith /*@ 9470b4b7b1cSBarry Smith SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` object to indicate the line search algorithm that should be used by a given `SNES` solver 948f40b411bSPeter Brune 949c3339decSBarry Smith Logically Collective 950f190f2fcSBarry Smith 951f40b411bSPeter Brune Input Parameters: 952420bcc1bSBarry Smith + linesearch - the line search context 953ceaaa498SBarry Smith - type - The type of line search to be used, see `SNESLineSearchType` 9541fe24845SBarry Smith 9553c7db156SBarry Smith Options Database Key: 956a99ef635SJonas Heinzmann . -snes_linesearch_type <type> - basic (or equivalently none), bt, secant, cp, nleqerr, bisection, shell 957cd7522eaSPeter Brune 958f40b411bSPeter Brune Level: intermediate 959f40b411bSPeter Brune 9600b4b7b1cSBarry Smith Note: 9610b4b7b1cSBarry Smith The `SNESLineSearch` object is generally obtained with `SNESGetLineSearch()` 9620b4b7b1cSBarry Smith 9630b4b7b1cSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`, 9640b4b7b1cSBarry Smith `SNESGetLineSearch()` 965f40b411bSPeter Brune @*/ 966d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 967d71ae5a4SJacob Faibussowitsch { 968bf7f4e0aSPeter Brune PetscBool match; 9695f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 970bf7f4e0aSPeter Brune 971bf7f4e0aSPeter Brune PetscFunctionBegin; 972f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9734f572ea9SToby Isaac PetscAssertPointer(type, 2); 974bf7f4e0aSPeter Brune 9759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 9763ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 977bf7f4e0aSPeter Brune 9789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 9796adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 980bf7f4e0aSPeter Brune /* Destroy the previous private line search context */ 981213acdd3SPierre Jolivet PetscTryTypeMethod(linesearch, destroy); 9820298fd71SBarry Smith linesearch->ops->destroy = NULL; 983f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9849e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9859e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9869e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9879e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 988bf7f4e0aSPeter Brune 9899566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 9909566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 9913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 992bf7f4e0aSPeter Brune } 993bf7f4e0aSPeter Brune 994f40b411bSPeter Brune /*@ 995f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 996f40b411bSPeter Brune 997f40b411bSPeter Brune Input Parameters: 998420bcc1bSBarry Smith + linesearch - the line search context 999420bcc1bSBarry Smith - snes - The `SNES` instance 1000f40b411bSPeter Brune 100178bcb3b5SPeter Brune Level: developer 100278bcb3b5SPeter Brune 1003f6dfbefdSBarry Smith Note: 1004f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 1005f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 100678bcb3b5SPeter Brune implementations. 1007f40b411bSPeter Brune 1008420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()` 1009f40b411bSPeter Brune @*/ 1010d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1011d71ae5a4SJacob Faibussowitsch { 1012bf7f4e0aSPeter Brune PetscFunctionBegin; 1013f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1014bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 1015bf7f4e0aSPeter Brune linesearch->snes = snes; 10163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1017bf7f4e0aSPeter Brune } 1018bf7f4e0aSPeter Brune 1019f40b411bSPeter Brune /*@ 1020f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 1021f6dfbefdSBarry Smith 1022f6dfbefdSBarry Smith Not Collective 1023f40b411bSPeter Brune 10242fe279fdSBarry Smith Input Parameter: 1025420bcc1bSBarry Smith . linesearch - the line search context 1026f40b411bSPeter Brune 10272fe279fdSBarry Smith Output Parameter: 10282fe279fdSBarry Smith . snes - The `SNES` instance 1029f40b411bSPeter Brune 10308141a3b9SPeter Brune Level: developer 1031f40b411bSPeter Brune 1032420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()` 1033f40b411bSPeter Brune @*/ 1034d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1035d71ae5a4SJacob Faibussowitsch { 1036bf7f4e0aSPeter Brune PetscFunctionBegin; 1037f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10384f572ea9SToby Isaac PetscAssertPointer(snes, 2); 1039bf7f4e0aSPeter Brune *snes = linesearch->snes; 10403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1041bf7f4e0aSPeter Brune } 1042bf7f4e0aSPeter Brune 1043f40b411bSPeter Brune /*@ 1044a99ef635SJonas Heinzmann SNESLineSearchGetLambda - Gets the last line search `lambda` used 1045f40b411bSPeter Brune 1046f6dfbefdSBarry Smith Not Collective 1047f6dfbefdSBarry Smith 10482fe279fdSBarry Smith Input Parameter: 1049420bcc1bSBarry Smith . linesearch - the line search context 1050f40b411bSPeter Brune 10512fe279fdSBarry Smith Output Parameter: 1052a99ef635SJonas Heinzmann . lambda - The last `lambda` (scaling of the solution udpate) computed during `SNESLineSearchApply()` 1053f40b411bSPeter Brune 105478bcb3b5SPeter Brune Level: advanced 105578bcb3b5SPeter Brune 1056f6dfbefdSBarry Smith Note: 10578e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 105878bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 1059f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the 1060a99ef635SJonas Heinzmann line search `lambda` using the argument -snes_qn_scaling ls. 106178bcb3b5SPeter Brune 1062420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1063f40b411bSPeter Brune @*/ 1064d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1065d71ae5a4SJacob Faibussowitsch { 10666a388c36SPeter Brune PetscFunctionBegin; 1067f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10684f572ea9SToby Isaac PetscAssertPointer(lambda, 2); 10696a388c36SPeter Brune *lambda = linesearch->lambda; 10703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10716a388c36SPeter Brune } 10726a388c36SPeter Brune 1073f40b411bSPeter Brune /*@ 1074a99ef635SJonas Heinzmann SNESLineSearchSetLambda - Sets the line search `lambda` (scaling of the solution update) 1075f40b411bSPeter Brune 1076f40b411bSPeter Brune Input Parameters: 10778e557f58SPeter Brune + linesearch - line search context 1078a99ef635SJonas Heinzmann - lambda - The `lambda` to use 1079f40b411bSPeter Brune 108020f4b53cSBarry Smith Level: advanced 108120f4b53cSBarry Smith 1082f6dfbefdSBarry Smith Note: 1083f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()` 1084a99ef635SJonas Heinzmann to set the final `lambda`. This routine (and `SNESLineSearchGetLambda()`) were 1085a99ef635SJonas Heinzmann added to facilitate Quasi-Newton methods that use the previous `lambda` 1086cd7522eaSPeter Brune as an inner scaling parameter. 1087cd7522eaSPeter Brune 1088420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()` 1089f40b411bSPeter Brune @*/ 1090d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1091d71ae5a4SJacob Faibussowitsch { 10926a388c36SPeter Brune PetscFunctionBegin; 1093f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10946a388c36SPeter Brune linesearch->lambda = lambda; 10953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10966a388c36SPeter Brune } 10976a388c36SPeter Brune 1098f40b411bSPeter Brune /*@ 1099ceaaa498SBarry Smith SNESLineSearchGetTolerances - Gets the tolerances for the line search. 1100f40b411bSPeter Brune 1101f6dfbefdSBarry Smith Not Collective 1102f6dfbefdSBarry Smith 1103f899ff85SJose E. Roman Input Parameter: 1104420bcc1bSBarry Smith . linesearch - the line search context 1105f40b411bSPeter Brune 1106f40b411bSPeter Brune Output Parameters: 1107a99ef635SJonas Heinzmann + minlambda - The minimum `lambda` allowed 1108a99ef635SJonas Heinzmann . maxlambda - The maximum `lambda` allowed 1109516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1110516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1111a99ef635SJonas Heinzmann . ltol - The change in `lambda` tolerance for iterative line searches 1112a99ef635SJonas Heinzmann - max_it - The maximum number of iterations of the line search 1113f40b411bSPeter Brune 111478bcb3b5SPeter Brune Level: intermediate 111578bcb3b5SPeter Brune 1116f6dfbefdSBarry Smith Note: 111778bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 11183c7d6663SPeter Brune the type requires. 1119516fe3c3SPeter Brune 1120420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1121f40b411bSPeter Brune @*/ 1122a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *minlambda, PetscReal *maxlambda, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_it) 1123d71ae5a4SJacob Faibussowitsch { 11246a388c36SPeter Brune PetscFunctionBegin; 1125f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1126a99ef635SJonas Heinzmann if (minlambda) { 1127a99ef635SJonas Heinzmann PetscAssertPointer(minlambda, 2); 1128a99ef635SJonas Heinzmann *minlambda = linesearch->minlambda; 1129516fe3c3SPeter Brune } 1130a99ef635SJonas Heinzmann if (maxlambda) { 1131a99ef635SJonas Heinzmann PetscAssertPointer(maxlambda, 3); 1132a99ef635SJonas Heinzmann *maxlambda = linesearch->maxlambda; 1133516fe3c3SPeter Brune } 1134516fe3c3SPeter Brune if (rtol) { 11354f572ea9SToby Isaac PetscAssertPointer(rtol, 4); 1136516fe3c3SPeter Brune *rtol = linesearch->rtol; 1137516fe3c3SPeter Brune } 1138516fe3c3SPeter Brune if (atol) { 11394f572ea9SToby Isaac PetscAssertPointer(atol, 5); 1140516fe3c3SPeter Brune *atol = linesearch->atol; 1141516fe3c3SPeter Brune } 1142516fe3c3SPeter Brune if (ltol) { 11434f572ea9SToby Isaac PetscAssertPointer(ltol, 6); 1144516fe3c3SPeter Brune *ltol = linesearch->ltol; 1145516fe3c3SPeter Brune } 1146a99ef635SJonas Heinzmann if (max_it) { 1147a99ef635SJonas Heinzmann PetscAssertPointer(max_it, 7); 1148a99ef635SJonas Heinzmann *max_it = linesearch->max_it; 1149516fe3c3SPeter Brune } 11503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11516a388c36SPeter Brune } 11526a388c36SPeter Brune 1153f40b411bSPeter Brune /*@ 1154ceaaa498SBarry Smith SNESLineSearchSetTolerances - Sets the tolerances for the linesearch. 1155f40b411bSPeter Brune 1156f6dfbefdSBarry Smith Collective 1157f6dfbefdSBarry Smith 1158f40b411bSPeter Brune Input Parameters: 1159420bcc1bSBarry Smith + linesearch - the line search context 1160a99ef635SJonas Heinzmann . minlambda - The minimum `lambda` allowed 1161a99ef635SJonas Heinzmann . maxlambda - The maximum `lambda` allowed 1162516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1163516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1164a99ef635SJonas Heinzmann . ltol - The change in `lambda` tolerance for iterative line searches 1165420bcc1bSBarry Smith - max_it - The maximum number of iterations of the line search 1166420bcc1bSBarry Smith 1167420bcc1bSBarry Smith Options Database Keys: 1168a99ef635SJonas Heinzmann + -snes_linesearch_minlambda - The minimum `lambda` allowed 1169a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda - The maximum `lambda` allowed 1170420bcc1bSBarry Smith . -snes_linesearch_rtol - Relative tolerance for iterative line searches 1171420bcc1bSBarry Smith . -snes_linesearch_atol - Absolute tolerance for iterative line searches 1172a99ef635SJonas Heinzmann . -snes_linesearch_ltol - Change in `lambda` tolerance for iterative line searches 1173420bcc1bSBarry Smith - -snes_linesearch_max_it - The number of iterations for iterative line searches 1174f40b411bSPeter Brune 117520f4b53cSBarry Smith Level: intermediate 117620f4b53cSBarry Smith 1177f6dfbefdSBarry Smith Note: 1178420bcc1bSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument. 1179f40b411bSPeter Brune 1180420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1181f40b411bSPeter Brune @*/ 1182a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal minlambda, PetscReal maxlambda, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it) 1183d71ae5a4SJacob Faibussowitsch { 11846a388c36SPeter Brune PetscFunctionBegin; 1185f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1186a99ef635SJonas Heinzmann PetscValidLogicalCollectiveReal(linesearch, minlambda, 2); 1187a99ef635SJonas Heinzmann PetscValidLogicalCollectiveReal(linesearch, maxlambda, 3); 1188d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1189d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1190d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1191420bcc1bSBarry Smith PetscValidLogicalCollectiveInt(linesearch, max_it, 7); 1192d3952184SSatish Balay 1193a99ef635SJonas Heinzmann if (minlambda != (PetscReal)PETSC_DEFAULT) { 1194a99ef635SJonas Heinzmann PetscCheck(minlambda >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum lambda %14.12e must be non-negative", (double)minlambda); 1195a99ef635SJonas Heinzmann PetscCheck(minlambda < maxlambda, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum lambda %14.12e must be smaller than maximum lambda %14.12e", (double)minlambda, (double)maxlambda); 1196a99ef635SJonas Heinzmann linesearch->minlambda = minlambda; 1197d3952184SSatish Balay } 1198d3952184SSatish Balay 1199a99ef635SJonas Heinzmann if (maxlambda != (PetscReal)PETSC_DEFAULT) { 1200a99ef635SJonas Heinzmann PetscCheck(maxlambda > 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum lambda %14.12e must be positive", (double)maxlambda); 1201a99ef635SJonas Heinzmann linesearch->maxlambda = maxlambda; 1202d3952184SSatish Balay } 1203d3952184SSatish Balay 120413bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) { 12052061ca32SJunchao 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); 1206516fe3c3SPeter Brune linesearch->rtol = rtol; 1207d3952184SSatish Balay } 1208d3952184SSatish Balay 120913bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) { 12105f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1211516fe3c3SPeter Brune linesearch->atol = atol; 1212d3952184SSatish Balay } 1213d3952184SSatish Balay 121413bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) { 12155f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1216516fe3c3SPeter Brune linesearch->ltol = ltol; 1217d3952184SSatish Balay } 1218d3952184SSatish Balay 1219420bcc1bSBarry Smith if (max_it != PETSC_DEFAULT) { 1220420bcc1bSBarry Smith PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it); 1221a99ef635SJonas Heinzmann linesearch->max_it = max_it; 1222d3952184SSatish Balay } 12233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12246a388c36SPeter Brune } 12256a388c36SPeter Brune 1226f40b411bSPeter Brune /*@ 1227f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1228f40b411bSPeter Brune 12292fe279fdSBarry Smith Input Parameter: 1230420bcc1bSBarry Smith . linesearch - the line search context 1231f40b411bSPeter Brune 12322fe279fdSBarry Smith Output Parameter: 12338e557f58SPeter Brune . damping - The damping parameter 1234f40b411bSPeter Brune 123578bcb3b5SPeter Brune Level: advanced 1236f40b411bSPeter Brune 1237420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN` 1238f40b411bSPeter Brune @*/ 1239d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1240d71ae5a4SJacob Faibussowitsch { 12416a388c36SPeter Brune PetscFunctionBegin; 1242f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12434f572ea9SToby Isaac PetscAssertPointer(damping, 2); 12446a388c36SPeter Brune *damping = linesearch->damping; 12453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12466a388c36SPeter Brune } 12476a388c36SPeter Brune 1248f40b411bSPeter Brune /*@ 1249fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1250f40b411bSPeter Brune 1251f40b411bSPeter Brune Input Parameters: 1252420bcc1bSBarry Smith + linesearch - the line search context 125303fd4120SBarry Smith - damping - The damping parameter 1254f40b411bSPeter Brune 1255f6dfbefdSBarry Smith Options Database Key: 1256f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value 1257f6dfbefdSBarry Smith 1258f40b411bSPeter Brune Level: intermediate 1259f40b411bSPeter Brune 1260f6dfbefdSBarry Smith Note: 1261f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1262a99ef635SJonas Heinzmann The use of the damping parameter in the `SNESLINESEARCHSECANT` and `SNESLINESEARCHCP` line searches is much more subtle; 1263a99ef635SJonas Heinzmann it is used as a starting point for the secant method. Depending on the choice for `maxlambda`, 1264a99ef635SJonas Heinzmann the eventual `lambda` may be greater than the damping parameter however. 1265a99ef635SJonas Heinzmann For `SNESLINESEARCHBISECTION` and `SNESLINESEARCHBT` the damping is instead used as the initial guess, 1266a99ef635SJonas Heinzmann below which the line search will not go. Hence, it is the maximum possible value for `lambda`. 1267cd7522eaSPeter Brune 1268420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()` 1269f40b411bSPeter Brune @*/ 1270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1271d71ae5a4SJacob Faibussowitsch { 12726a388c36SPeter Brune PetscFunctionBegin; 1273f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12746a388c36SPeter Brune linesearch->damping = damping; 12753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12766a388c36SPeter Brune } 12776a388c36SPeter Brune 127859405d5eSPeter Brune /*@ 127959405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 128059405d5eSPeter Brune 1281f6dfbefdSBarry Smith Input Parameter: 1282420bcc1bSBarry Smith . linesearch - the line search context 128359405d5eSPeter Brune 1284f6dfbefdSBarry Smith Output Parameter: 12858e557f58SPeter Brune . order - The order 128659405d5eSPeter Brune 128759405d5eSPeter Brune Level: intermediate 128859405d5eSPeter Brune 1289420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()` 129059405d5eSPeter Brune @*/ 1291d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1292d71ae5a4SJacob Faibussowitsch { 129359405d5eSPeter Brune PetscFunctionBegin; 129459405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12954f572ea9SToby Isaac PetscAssertPointer(order, 2); 129659405d5eSPeter Brune *order = linesearch->order; 12973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129859405d5eSPeter Brune } 129959405d5eSPeter Brune 130059405d5eSPeter Brune /*@ 13011f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 130259405d5eSPeter Brune 130359405d5eSPeter Brune Input Parameters: 1304420bcc1bSBarry Smith + linesearch - the line search context 1305ceaaa498SBarry Smith - order - The order 130659405d5eSPeter Brune 130759405d5eSPeter Brune Level: intermediate 130859405d5eSPeter Brune 1309ceaaa498SBarry Smith Values for `order`\: 1310f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1311f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1312f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 131378bcb3b5SPeter Brune 1314420bcc1bSBarry Smith Options Database Key: 1315420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3) 1316420bcc1bSBarry Smith 1317ceaaa498SBarry Smith Note: 1318ceaaa498SBarry Smith These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP` 131959405d5eSPeter Brune 1320420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 132159405d5eSPeter Brune @*/ 1322d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1323d71ae5a4SJacob Faibussowitsch { 132459405d5eSPeter Brune PetscFunctionBegin; 132559405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 132659405d5eSPeter Brune linesearch->order = order; 13273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132859405d5eSPeter Brune } 132959405d5eSPeter Brune 1330f40b411bSPeter Brune /*@ 1331420bcc1bSBarry Smith SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1332f40b411bSPeter Brune 1333f6dfbefdSBarry Smith Not Collective 1334f6dfbefdSBarry Smith 1335f899ff85SJose E. Roman Input Parameter: 1336420bcc1bSBarry Smith . linesearch - the line search context 1337f40b411bSPeter Brune 1338f40b411bSPeter Brune Output Parameters: 1339f40b411bSPeter Brune + xnorm - The norm of the current solution 1340a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1341a99ef635SJonas Heinzmann - ynorm - The norm of the current update (after scaling by the linesearch computed `lambda`) 1342f40b411bSPeter Brune 134378bcb3b5SPeter Brune Level: developer 1344f40b411bSPeter Brune 1345420bcc1bSBarry Smith Notes: 1346420bcc1bSBarry Smith Some values may not be up-to-date at particular points in the code. 1347a68bbae5SBarry Smith 1348a68bbae5SBarry Smith This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1349a68bbae5SBarry Smith computed values. 1350a68bbae5SBarry Smith 13510241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1352f40b411bSPeter Brune @*/ 1353d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1354d71ae5a4SJacob Faibussowitsch { 1355bf7f4e0aSPeter Brune PetscFunctionBegin; 1356f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1357f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1358f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1359f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 13603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1361bf7f4e0aSPeter Brune } 1362bf7f4e0aSPeter Brune 1363f40b411bSPeter Brune /*@ 13641bd63e3eSSatish Balay SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1365f40b411bSPeter Brune 1366c3339decSBarry Smith Collective 1367f6dfbefdSBarry Smith 1368f40b411bSPeter Brune Input Parameters: 1369420bcc1bSBarry Smith + linesearch - the line search context 1370f40b411bSPeter Brune . xnorm - The norm of the current solution 1371a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1372a99ef635SJonas Heinzmann - ynorm - The norm of the current update (after scaling by the linesearch computed `lambda`) 1373f40b411bSPeter Brune 1374f6dfbefdSBarry Smith Level: developer 1375f40b411bSPeter Brune 1376420bcc1bSBarry Smith Note: 1377420bcc1bSBarry Smith This is called by the line search routines to store the values they have just computed 1378420bcc1bSBarry Smith 1379420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1380f40b411bSPeter Brune @*/ 1381d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1382d71ae5a4SJacob Faibussowitsch { 13836a388c36SPeter Brune PetscFunctionBegin; 1384f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 13856a388c36SPeter Brune linesearch->xnorm = xnorm; 13866a388c36SPeter Brune linesearch->fnorm = fnorm; 13876a388c36SPeter Brune linesearch->ynorm = ynorm; 13883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13896a388c36SPeter Brune } 13906a388c36SPeter Brune 1391f40b411bSPeter Brune /*@ 1392420bcc1bSBarry Smith SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`. 1393f40b411bSPeter Brune 13942fe279fdSBarry Smith Input Parameter: 1395420bcc1bSBarry Smith . linesearch - the line search context 1396f40b411bSPeter Brune 139720f4b53cSBarry Smith Options Database Key: 13988e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1399f40b411bSPeter Brune 1400f40b411bSPeter Brune Level: intermediate 1401f40b411bSPeter Brune 1402420bcc1bSBarry Smith Developer Note: 1403420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms 1404420bcc1bSBarry Smith 1405420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1406f40b411bSPeter Brune @*/ 1407d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1408d71ae5a4SJacob Faibussowitsch { 14099bd66eb0SPeter Brune SNES snes; 1410bf388a1fSBarry Smith 14116a388c36SPeter Brune PetscFunctionBegin; 14126a388c36SPeter Brune if (linesearch->norms) { 14139bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 14149566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 14159566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14169566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14179566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 14189bd66eb0SPeter Brune } else { 14199566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 14209566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14219566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14229566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 14239566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 14249566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 14256a388c36SPeter Brune } 14269bd66eb0SPeter Brune } 14273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14286a388c36SPeter Brune } 14296a388c36SPeter Brune 14306f263ca3SPeter Brune /*@ 14316f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 14326f263ca3SPeter Brune 14336f263ca3SPeter Brune Input Parameters: 1434420bcc1bSBarry Smith + linesearch - the line search context 143578bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 14366f263ca3SPeter Brune 143720f4b53cSBarry Smith Options Database Key: 1438420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search 14396f263ca3SPeter Brune 144020f4b53cSBarry Smith Level: intermediate 144120f4b53cSBarry Smith 1442f6dfbefdSBarry Smith Note: 1443f6dfbefdSBarry 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. 14446f263ca3SPeter Brune 1445420bcc1bSBarry Smith Developer Note: 1446420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms 1447420bcc1bSBarry Smith 1448420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 14496f263ca3SPeter Brune @*/ 1450d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1451d71ae5a4SJacob Faibussowitsch { 14526f263ca3SPeter Brune PetscFunctionBegin; 14536f263ca3SPeter Brune linesearch->norms = flg; 14543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14556f263ca3SPeter Brune } 14566f263ca3SPeter Brune 1457f40b411bSPeter Brune /*@ 1458f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1459f6dfbefdSBarry Smith 14608f14a041SBarry Smith Not Collective but the vectors are parallel 1461f40b411bSPeter Brune 1462f899ff85SJose E. Roman Input Parameter: 1463420bcc1bSBarry Smith . linesearch - the line search context 1464f40b411bSPeter Brune 1465f40b411bSPeter Brune Output Parameters: 14666232e825SPeter Brune + X - Solution vector 14676232e825SPeter Brune . F - Function vector 14686232e825SPeter Brune . Y - Search direction vector 14696232e825SPeter Brune . W - Solution work vector 14706232e825SPeter Brune - G - Function work vector 14716232e825SPeter Brune 147220f4b53cSBarry Smith Level: advanced 147320f4b53cSBarry Smith 14747bba9028SPeter Brune Notes: 147520f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a 147620f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the 147720f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the 14786232e825SPeter Brune function evaluated at the new solution. 1479f40b411bSPeter Brune 1480f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 14812a7a6963SBarry Smith 1482420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1483f40b411bSPeter Brune @*/ 1484d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1485d71ae5a4SJacob Faibussowitsch { 14866a388c36SPeter Brune PetscFunctionBegin; 1487f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14886a388c36SPeter Brune if (X) { 14894f572ea9SToby Isaac PetscAssertPointer(X, 2); 14906a388c36SPeter Brune *X = linesearch->vec_sol; 14916a388c36SPeter Brune } 14926a388c36SPeter Brune if (F) { 14934f572ea9SToby Isaac PetscAssertPointer(F, 3); 14946a388c36SPeter Brune *F = linesearch->vec_func; 14956a388c36SPeter Brune } 14966a388c36SPeter Brune if (Y) { 14974f572ea9SToby Isaac PetscAssertPointer(Y, 4); 14986a388c36SPeter Brune *Y = linesearch->vec_update; 14996a388c36SPeter Brune } 15006a388c36SPeter Brune if (W) { 15014f572ea9SToby Isaac PetscAssertPointer(W, 5); 15026a388c36SPeter Brune *W = linesearch->vec_sol_new; 15036a388c36SPeter Brune } 15046a388c36SPeter Brune if (G) { 15054f572ea9SToby Isaac PetscAssertPointer(G, 6); 15066a388c36SPeter Brune *G = linesearch->vec_func_new; 15076a388c36SPeter Brune } 15083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15096a388c36SPeter Brune } 15106a388c36SPeter Brune 1511f40b411bSPeter Brune /*@ 1512f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1513f6dfbefdSBarry Smith 1514c3339decSBarry Smith Logically Collective 1515f40b411bSPeter Brune 1516f40b411bSPeter Brune Input Parameters: 1517420bcc1bSBarry Smith + linesearch - the line search context 15186232e825SPeter Brune . X - Solution vector 15196232e825SPeter Brune . F - Function vector 15206232e825SPeter Brune . Y - Search direction vector 15216232e825SPeter Brune . W - Solution work vector 15226232e825SPeter Brune - G - Function work vector 1523f40b411bSPeter Brune 1524420bcc1bSBarry Smith Level: developer 1525f40b411bSPeter Brune 1526420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1527f40b411bSPeter Brune @*/ 1528d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1529d71ae5a4SJacob Faibussowitsch { 15306a388c36SPeter Brune PetscFunctionBegin; 1531f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15326a388c36SPeter Brune if (X) { 15336a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 15346a388c36SPeter Brune linesearch->vec_sol = X; 15356a388c36SPeter Brune } 15366a388c36SPeter Brune if (F) { 15376a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 15386a388c36SPeter Brune linesearch->vec_func = F; 15396a388c36SPeter Brune } 15406a388c36SPeter Brune if (Y) { 15416a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 15426a388c36SPeter Brune linesearch->vec_update = Y; 15436a388c36SPeter Brune } 15446a388c36SPeter Brune if (W) { 15456a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 15466a388c36SPeter Brune linesearch->vec_sol_new = W; 15476a388c36SPeter Brune } 15486a388c36SPeter Brune if (G) { 15496a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 15506a388c36SPeter Brune linesearch->vec_func_new = G; 15516a388c36SPeter Brune } 15523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15536a388c36SPeter Brune } 15546a388c36SPeter Brune 1555cc4c1da9SBarry Smith /*@ 1556f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1557f6dfbefdSBarry Smith `SNESLineSearch` options in the database. 1558e7058c64SPeter Brune 1559c3339decSBarry Smith Logically Collective 1560e7058c64SPeter Brune 1561e7058c64SPeter Brune Input Parameters: 1562f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1563e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1564e7058c64SPeter Brune 156520f4b53cSBarry Smith Level: advanced 156620f4b53cSBarry Smith 1567f6dfbefdSBarry Smith Note: 1568e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1569e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1570e7058c64SPeter Brune 1571420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1572e7058c64SPeter Brune @*/ 1573d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1574d71ae5a4SJacob Faibussowitsch { 1575e7058c64SPeter Brune PetscFunctionBegin; 1576f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15779566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 15783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1579e7058c64SPeter Brune } 1580e7058c64SPeter Brune 1581cc4c1da9SBarry Smith /*@ 1582f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1583f1c6b773SPeter Brune SNESLineSearch options in the database. 1584e7058c64SPeter Brune 1585e7058c64SPeter Brune Not Collective 1586e7058c64SPeter Brune 1587e7058c64SPeter Brune Input Parameter: 1588f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 1589e7058c64SPeter Brune 1590e7058c64SPeter Brune Output Parameter: 1591e7058c64SPeter Brune . prefix - pointer to the prefix string used 1592e7058c64SPeter Brune 1593e7058c64SPeter Brune Level: advanced 1594e7058c64SPeter Brune 1595420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1596e7058c64SPeter Brune @*/ 1597d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1598d71ae5a4SJacob Faibussowitsch { 1599e7058c64SPeter Brune PetscFunctionBegin; 1600f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 16023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1603e7058c64SPeter Brune } 1604bf7f4e0aSPeter Brune 16058d359177SBarry Smith /*@C 1606f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1607f40b411bSPeter Brune 1608d8d19677SJose E. Roman Input Parameters: 1609f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1610f40b411bSPeter Brune - nwork - the number of work vectors 1611f40b411bSPeter Brune 1612f40b411bSPeter Brune Level: developer 1613f40b411bSPeter Brune 1614420bcc1bSBarry Smith Developer Note: 1615420bcc1bSBarry Smith This is called from within the set up routines for each of the line search types `SNESLineSearchType` 1616420bcc1bSBarry Smith 1617420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()` 1618f40b411bSPeter Brune @*/ 1619d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1620d71ae5a4SJacob Faibussowitsch { 1621bf7f4e0aSPeter Brune PetscFunctionBegin; 16220fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 16239566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 16243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1625bf7f4e0aSPeter Brune } 1626bf7f4e0aSPeter Brune 1627f40b411bSPeter Brune /*@ 1628422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1629f40b411bSPeter Brune 16302fe279fdSBarry Smith Input Parameter: 1631420bcc1bSBarry Smith . linesearch - the line search context 1632f40b411bSPeter Brune 16332fe279fdSBarry Smith Output Parameter: 1634422a814eSBarry Smith . result - The success or failure status 1635f40b411bSPeter Brune 163620f4b53cSBarry Smith Level: developer 163720f4b53cSBarry Smith 1638f6dfbefdSBarry Smith Note: 1639420bcc1bSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed 1640420bcc1bSBarry Smith (and set into the `SNES` convergence accordingly). 1641cd7522eaSPeter Brune 1642420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1643f40b411bSPeter Brune @*/ 1644d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1645d71ae5a4SJacob Faibussowitsch { 1646bf7f4e0aSPeter Brune PetscFunctionBegin; 1647f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16484f572ea9SToby Isaac PetscAssertPointer(result, 2); 1649422a814eSBarry Smith *result = linesearch->result; 16503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1651bf7f4e0aSPeter Brune } 1652bf7f4e0aSPeter Brune 1653f40b411bSPeter Brune /*@ 1654420bcc1bSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the line search application 1655f40b411bSPeter Brune 16565d83a8b1SBarry Smith Logically Collective; No Fortran Support 16575d83a8b1SBarry Smith 1658f40b411bSPeter Brune Input Parameters: 1659420bcc1bSBarry Smith + linesearch - the line search context 1660422a814eSBarry Smith - result - The success or failure status 1661f40b411bSPeter Brune 166220f4b53cSBarry Smith Level: developer 166320f4b53cSBarry Smith 1664f6dfbefdSBarry Smith Note: 1665420bcc1bSBarry Smith This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1666cd7522eaSPeter Brune the success or failure of the line search method. 1667cd7522eaSPeter Brune 1668420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1669f40b411bSPeter Brune @*/ 1670d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1671d71ae5a4SJacob Faibussowitsch { 16726a388c36SPeter Brune PetscFunctionBegin; 16735fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1674422a814eSBarry Smith linesearch->result = result; 16753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16766a388c36SPeter Brune } 16776a388c36SPeter Brune 1678ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 16799bd66eb0SPeter Brune /*@C 1680f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16819bd66eb0SPeter Brune 1682c3339decSBarry Smith Logically Collective 1683f6dfbefdSBarry Smith 16849bd66eb0SPeter Brune Input Parameters: 1685ceaaa498SBarry Smith + linesearch - the linesearch object 16868434afd1SBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence 1687d5def619SJonas Heinzmann . normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence 1688d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence 16899bd66eb0SPeter Brune 1690f6dfbefdSBarry Smith Level: advanced 16919bd66eb0SPeter Brune 1692ceaaa498SBarry Smith Notes: 169320f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 169420f4b53cSBarry Smith 169520f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated 169620f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`. 169720f4b53cSBarry Smith 1698d5def619SJonas Heinzmann The VI solvers further require special evaluation of the directional derivative (when assuming that there exists some $G(x)$ 1699d5def619SJonas Heinzmann for which the `SNESFunctionFn` $F(x) = grad G(x)$) such that it is only calculated on the inactive set. 1700d5def619SJonas Heinzmann This should be implemented by `dirderivfunc`. 1701d5def619SJonas Heinzmann 17029bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, 1703d5def619SJonas Heinzmann `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`, `SNESLineSearchVIDirDerivFn` 17049bd66eb0SPeter Brune @*/ 1705d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc, SNESLineSearchVIDirDerivFn *dirderivfunc) 1706d71ae5a4SJacob Faibussowitsch { 17079bd66eb0SPeter Brune PetscFunctionBegin; 1708f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 17099bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 17109bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 1711d5def619SJonas Heinzmann if (dirderivfunc) linesearch->ops->vidirderiv = dirderivfunc; 17123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17139bd66eb0SPeter Brune } 17149bd66eb0SPeter Brune 17159bd66eb0SPeter Brune /*@C 1716f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 17179bd66eb0SPeter Brune 1718f6dfbefdSBarry Smith Not Collective 1719f6dfbefdSBarry Smith 1720f899ff85SJose E. Roman Input Parameter: 1721f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()` 17229bd66eb0SPeter Brune 17239bd66eb0SPeter Brune Output Parameters: 17248434afd1SBarry Smith + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence 1725d5def619SJonas Heinzmann . normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence 1726d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence 17279bd66eb0SPeter Brune 1728f6dfbefdSBarry Smith Level: advanced 17299bd66eb0SPeter Brune 17309bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 17318434afd1SBarry Smith `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn` 17329bd66eb0SPeter Brune @*/ 1733d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc, SNESLineSearchVIDirDerivFn **dirderivfunc) 1734d71ae5a4SJacob Faibussowitsch { 17359bd66eb0SPeter Brune PetscFunctionBegin; 17369bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 17379bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 1738d5def619SJonas Heinzmann if (dirderivfunc) *dirderivfunc = linesearch->ops->vidirderiv; 17393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17409bd66eb0SPeter Brune } 17419bd66eb0SPeter Brune 1742bf7f4e0aSPeter Brune /*@C 1743420bcc1bSBarry Smith SNESLineSearchRegister - register a line search type `SNESLineSearchType` 1744ceaaa498SBarry Smith 1745cc4c1da9SBarry Smith Logically Collective, No Fortran Support 1746cc4c1da9SBarry Smith 1747ceaaa498SBarry Smith Input Parameters: 1748420bcc1bSBarry Smith + sname - name of the `SNESLineSearchType()` 1749ceaaa498SBarry Smith - function - the creation function for that type 1750ceaaa498SBarry Smith 1751ceaaa498SBarry Smith Calling sequence of `function`: 1752ceaaa498SBarry Smith . ls - the line search context 1753bf7f4e0aSPeter Brune 1754bf7f4e0aSPeter Brune Level: advanced 1755f6dfbefdSBarry Smith 1756420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1757bf7f4e0aSPeter Brune @*/ 1758ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls)) 1759d71ae5a4SJacob Faibussowitsch { 1760bf7f4e0aSPeter Brune PetscFunctionBegin; 17619566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 17629566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 17633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1764bf7f4e0aSPeter Brune } 1765