1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2bf7f4e0aSPeter Brune 3f1c6b773SPeter Brune PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 40298fd71SBarry Smith PetscFunctionList SNESLineSearchList = NULL; 5bf7f4e0aSPeter Brune 6f1c6b773SPeter Brune PetscClassId SNESLINESEARCH_CLASSID; 757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply; 8bf7f4e0aSPeter Brune 9dcf2fd19SBarry Smith /*@ 10f6dfbefdSBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object. 11dcf2fd19SBarry Smith 12c3339decSBarry Smith Logically Collective 13dcf2fd19SBarry Smith 142fe279fdSBarry Smith Input Parameter: 15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context 16dcf2fd19SBarry Smith 17dcf2fd19SBarry Smith Options Database Key: 18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 19f6dfbefdSBarry Smith into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those 20dcf2fd19SBarry Smith set via the options database 21dcf2fd19SBarry Smith 2220f4b53cSBarry Smith Level: advanced 2320f4b53cSBarry Smith 24dcf2fd19SBarry Smith Notes: 25f6dfbefdSBarry Smith There is no way to clear one specific monitor from a `SNESLineSearch` object. 26dcf2fd19SBarry Smith 27f6dfbefdSBarry Smith This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel 28dcf2fd19SBarry Smith that one. 29dcf2fd19SBarry Smith 30f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()` 31dcf2fd19SBarry Smith @*/ 32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 33d71ae5a4SJacob Faibussowitsch { 34dcf2fd19SBarry Smith PetscInt i; 35dcf2fd19SBarry Smith 36dcf2fd19SBarry Smith PetscFunctionBegin; 37dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 38dcf2fd19SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 3948a46eb9SPierre Jolivet if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i])); 40dcf2fd19SBarry Smith } 41dcf2fd19SBarry Smith ls->numbermonitors = 0; 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43dcf2fd19SBarry Smith } 44dcf2fd19SBarry Smith 45dcf2fd19SBarry Smith /*@ 46dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 47dcf2fd19SBarry Smith 48c3339decSBarry Smith Collective 49dcf2fd19SBarry Smith 502fe279fdSBarry Smith Input Parameter: 51dcf2fd19SBarry Smith . ls - the linesearch object 52dcf2fd19SBarry Smith 532fe279fdSBarry Smith Level: developer 542fe279fdSBarry Smith 55f6dfbefdSBarry Smith Note: 5620f4b53cSBarry Smith This routine is called by the `SNES` implementations. 57dcf2fd19SBarry Smith It does not typically need to be called by the user. 58dcf2fd19SBarry Smith 59f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()` 60dcf2fd19SBarry Smith @*/ 61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 62d71ae5a4SJacob Faibussowitsch { 63dcf2fd19SBarry Smith PetscInt i, n = ls->numbermonitors; 64dcf2fd19SBarry Smith 65dcf2fd19SBarry Smith PetscFunctionBegin; 6648a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i])); 673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68dcf2fd19SBarry Smith } 69dcf2fd19SBarry Smith 70dcf2fd19SBarry Smith /*@C 71dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 72dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 73dcf2fd19SBarry Smith progress. 74dcf2fd19SBarry Smith 75c3339decSBarry Smith Logically Collective 76dcf2fd19SBarry Smith 77dcf2fd19SBarry Smith Input Parameters: 78f6dfbefdSBarry Smith + ls - the `SNESLineSearch` context 79dcf2fd19SBarry Smith . f - the monitor function 80dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 8120f4b53cSBarry Smith monitor routine (use `NULL` if no context is desired) 82dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 8320f4b53cSBarry Smith (may be `NULL`) 8420f4b53cSBarry Smith 8520f4b53cSBarry Smith Level: intermediate 86dcf2fd19SBarry Smith 87f6dfbefdSBarry Smith Note: 88dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 89f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` multiple times; all will be called in the 90dcf2fd19SBarry Smith order in which they were set. 91dcf2fd19SBarry Smith 92e4094ef1SJacob Faibussowitsch Fortran Notes: 93f6dfbefdSBarry Smith Only a single monitor function can be set for each `SNESLineSearch` object 94dcf2fd19SBarry Smith 95f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()` 96dcf2fd19SBarry Smith @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **)) 98d71ae5a4SJacob Faibussowitsch { 9978064530SBarry Smith PetscInt i; 10078064530SBarry Smith PetscBool identical; 10178064530SBarry Smith 102dcf2fd19SBarry Smith PetscFunctionBegin; 103dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 10478064530SBarry Smith for (i = 0; i < ls->numbermonitors; i++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical)); 1063ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 10778064530SBarry Smith } 1085f80ce2aSJacob Faibussowitsch PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 109dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 110dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 111dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void *)mctx; 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113dcf2fd19SBarry Smith } 114dcf2fd19SBarry Smith 115dcf2fd19SBarry Smith /*@C 116f6dfbefdSBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries 117dcf2fd19SBarry Smith 118c3339decSBarry Smith Collective 119dcf2fd19SBarry Smith 120dcf2fd19SBarry Smith Input Parameters: 121f6dfbefdSBarry Smith + ls - the `SNES` linesearch object 122f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat` 123dcf2fd19SBarry Smith 124dcf2fd19SBarry Smith Level: intermediate 125dcf2fd19SBarry Smith 126f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()` 127dcf2fd19SBarry Smith @*/ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) 129d71ae5a4SJacob Faibussowitsch { 130d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 131dcf2fd19SBarry Smith Vec Y, W, G; 132dcf2fd19SBarry Smith 133dcf2fd19SBarry Smith PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n")); 1379566063dSJacob Faibussowitsch PetscCall(VecView(Y, viewer)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n")); 1399566063dSJacob Faibussowitsch PetscCall(VecView(W, viewer)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n")); 1419566063dSJacob Faibussowitsch PetscCall(VecView(G, viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144dcf2fd19SBarry Smith } 145dcf2fd19SBarry Smith 146f40b411bSPeter Brune /*@ 147cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 148f40b411bSPeter Brune 149f6dfbefdSBarry Smith Logically Collective 150f40b411bSPeter Brune 1512fe279fdSBarry Smith Input Parameter: 152f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context). 153f40b411bSPeter Brune 1542fe279fdSBarry Smith Output Parameter: 1558e557f58SPeter Brune . outlinesearch - the new linesearch context 156f40b411bSPeter Brune 157162e0bf5SPeter Brune Level: developer 158162e0bf5SPeter Brune 159f6dfbefdSBarry Smith Note: 160f6dfbefdSBarry Smith The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance 161f6dfbefdSBarry Smith already associated with the `SNES`. 162f40b411bSPeter Brune 163f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()` 164f40b411bSPeter Brune @*/ 165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 166d71ae5a4SJacob Faibussowitsch { 167f1c6b773SPeter Brune SNESLineSearch linesearch; 168bf388a1fSBarry Smith 169bf7f4e0aSPeter Brune PetscFunctionBegin; 1704f572ea9SToby Isaac PetscAssertPointer(outlinesearch, 2); 1719566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 1720298fd71SBarry Smith *outlinesearch = NULL; 173f5af7f23SKarl Rupp 1749566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView)); 175bf7f4e0aSPeter Brune 1760298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1770298fd71SBarry Smith linesearch->vec_func_new = NULL; 1780298fd71SBarry Smith linesearch->vec_sol = NULL; 1790298fd71SBarry Smith linesearch->vec_func = NULL; 1800298fd71SBarry Smith linesearch->vec_update = NULL; 1819bd66eb0SPeter Brune 182bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 183bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 184bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 185bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 186422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 187bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 188bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 189bf7f4e0aSPeter Brune linesearch->damping = 1.0; 190bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 191bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 192516fe3c3SPeter Brune linesearch->rtol = 1e-8; 193516fe3c3SPeter Brune linesearch->atol = 1e-15; 194516fe3c3SPeter Brune linesearch->ltol = 1e-8; 1950298fd71SBarry Smith linesearch->precheckctx = NULL; 1960298fd71SBarry Smith linesearch->postcheckctx = NULL; 19759405d5eSPeter Brune linesearch->max_its = 1; 198bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 1993add74b1SBarry Smith linesearch->monitor = NULL; 200bf7f4e0aSPeter Brune *outlinesearch = linesearch; 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 202bf7f4e0aSPeter Brune } 203bf7f4e0aSPeter Brune 204f40b411bSPeter Brune /*@ 20578bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 20678bcb3b5SPeter Brune any required vectors. 207f40b411bSPeter Brune 208c3339decSBarry Smith Collective 209f40b411bSPeter Brune 2102fe279fdSBarry Smith Input Parameter: 211f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 212f40b411bSPeter Brune 21320f4b53cSBarry Smith Level: advanced 21420f4b53cSBarry Smith 215f6dfbefdSBarry Smith Note: 216f6dfbefdSBarry Smith For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`. 217cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 21878bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 219f6dfbefdSBarry Smith of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be 220cd7522eaSPeter Brune allocated upfront. 221cd7522eaSPeter Brune 222f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()` 223f40b411bSPeter Brune @*/ 224d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 225d71ae5a4SJacob Faibussowitsch { 226bf7f4e0aSPeter Brune PetscFunctionBegin; 22748a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 228bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 22948a46eb9SPierre Jolivet if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 23048a46eb9SPierre Jolivet if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 231dbbe0bcdSBarry Smith if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup); 2329566063dSJacob Faibussowitsch if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction)); 233bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 234bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 235bf7f4e0aSPeter Brune } 2363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 237bf7f4e0aSPeter Brune } 238bf7f4e0aSPeter Brune 239f40b411bSPeter Brune /*@ 240f6dfbefdSBarry Smith SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search. 241f40b411bSPeter Brune 242c3339decSBarry Smith Collective 243f40b411bSPeter Brune 2442fe279fdSBarry Smith Input Parameter: 245f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance. 246f40b411bSPeter Brune 24720f4b53cSBarry Smith Level: developer 24820f4b53cSBarry Smith 249f6dfbefdSBarry Smith Note: 250f6dfbefdSBarry Smith Usually only called by `SNESReset()` 251f190f2fcSBarry Smith 252f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 253f40b411bSPeter Brune @*/ 254d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 255d71ae5a4SJacob Faibussowitsch { 256bf7f4e0aSPeter Brune PetscFunctionBegin; 257dbbe0bcdSBarry Smith if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset); 258f5af7f23SKarl Rupp 2599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_sol_new)); 2609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_func_new)); 261bf7f4e0aSPeter Brune 2629566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 263f5af7f23SKarl Rupp 264bf7f4e0aSPeter Brune linesearch->nwork = 0; 265bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 267bf7f4e0aSPeter Brune } 268bf7f4e0aSPeter Brune 269ed07d7d7SPeter Brune /*@C 270f6dfbefdSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search 271f6dfbefdSBarry Smith ` 272e4094ef1SJacob Faibussowitsch 273ed07d7d7SPeter Brune Input Parameters: 274e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context 275e4094ef1SJacob Faibussowitsch - func - function evaluation routine, this is usually the function provided with `SNESSetFunction()` 276ed07d7d7SPeter Brune 277ed07d7d7SPeter Brune Level: developer 278ed07d7d7SPeter Brune 279f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()` 280ed07d7d7SPeter Brune @*/ 281d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec)) 282d71ae5a4SJacob Faibussowitsch { 283ed07d7d7SPeter Brune PetscFunctionBegin; 284ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 285ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 287ed07d7d7SPeter Brune } 288ed07d7d7SPeter Brune 28986d74e61SPeter Brune /*@C 290f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 291df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 292f190f2fcSBarry Smith determined the search direction. 29386d74e61SPeter Brune 294c3339decSBarry Smith Logically Collective 29586d74e61SPeter Brune 29686d74e61SPeter Brune Input Parameters: 297f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 29820f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPreCheck()` 29920f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 30086d74e61SPeter Brune 30186d74e61SPeter Brune Level: intermediate 30286d74e61SPeter Brune 303f6dfbefdSBarry Smith Note: 304f0b84518SBarry Smith Use `SNESLineSearchSetPostCheck()` to change the step after the line search. 305f0b84518SBarry Smith search is complete. 306f0b84518SBarry Smith 307f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 308f0b84518SBarry Smith 309f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 310869ba2dcSBarry Smith `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 311f0b84518SBarry Smith 31286d74e61SPeter Brune @*/ 313d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx) 314d71ae5a4SJacob Faibussowitsch { 3159bd66eb0SPeter Brune PetscFunctionBegin; 316f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 317f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 31886d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32086d74e61SPeter Brune } 32186d74e61SPeter Brune 32286d74e61SPeter Brune /*@C 32353deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 32486d74e61SPeter Brune 325f899ff85SJose E. Roman Input Parameter: 326f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 32786d74e61SPeter Brune 32886d74e61SPeter Brune Output Parameters: 32920f4b53cSBarry Smith + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchPreCheck` 33020f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 33186d74e61SPeter Brune 33286d74e61SPeter Brune Level: intermediate 33386d74e61SPeter Brune 334e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 33586d74e61SPeter Brune @*/ 336d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) 337d71ae5a4SJacob Faibussowitsch { 3389bd66eb0SPeter Brune PetscFunctionBegin; 339f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 340f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 34186d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34386d74e61SPeter Brune } 34486d74e61SPeter Brune 34586d74e61SPeter Brune /*@C 346f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 347f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 34886d74e61SPeter Brune 34920f4b53cSBarry Smith Logically Collective 35086d74e61SPeter Brune 35186d74e61SPeter Brune Input Parameters: 352f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 35320f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPostCheck` 35420f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 35586d74e61SPeter Brune 35686d74e61SPeter Brune Level: intermediate 35786d74e61SPeter Brune 358f0b84518SBarry Smith Notes: 359f0b84518SBarry Smith Use `SNESLineSearchSetPreCheck()` to change the step before the line search. 360f0b84518SBarry Smith search is complete. 361f0b84518SBarry Smith 362f0b84518SBarry Smith Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 363f0b84518SBarry Smith 364f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`, 365e4094ef1SJacob Faibussowitsch `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 36686d74e61SPeter Brune @*/ 367d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 368d71ae5a4SJacob Faibussowitsch { 36986d74e61SPeter Brune PetscFunctionBegin; 370f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 371f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 37286d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 3733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37486d74e61SPeter Brune } 37586d74e61SPeter Brune 37686d74e61SPeter Brune /*@C 377f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 37886d74e61SPeter Brune 379f899ff85SJose E. Roman Input Parameter: 380f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 38186d74e61SPeter Brune 38286d74e61SPeter Brune Output Parameters: 383f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck` 38420f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 38586d74e61SPeter Brune 38686d74e61SPeter Brune Level: intermediate 38786d74e61SPeter Brune 388f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 38986d74e61SPeter Brune @*/ 390d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 391d71ae5a4SJacob Faibussowitsch { 3929bd66eb0SPeter Brune PetscFunctionBegin; 393f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 394f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 39586d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39786d74e61SPeter Brune } 39886d74e61SPeter Brune 399f40b411bSPeter Brune /*@ 400f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 401f40b411bSPeter Brune 402c3339decSBarry Smith Logically Collective 403f40b411bSPeter Brune 404f40b411bSPeter Brune Input Parameters: 4057b1df9c1SPeter Brune + linesearch - The linesearch instance. 4067b1df9c1SPeter Brune . X - The current solution 4077b1df9c1SPeter Brune - Y - The step direction 408f40b411bSPeter Brune 4092fe279fdSBarry Smith Output Parameter: 4108e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 411f40b411bSPeter Brune 412f0b84518SBarry Smith Level: advanced 413f40b411bSPeter Brune 414f6dfbefdSBarry Smith Note: 415f6dfbefdSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` 416f6dfbefdSBarry Smith 417f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, 418*fe8e7dddSPierre Jolivet `SNESLineSearchGetPostCheck()` 419f40b411bSPeter Brune @*/ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) 421d71ae5a4SJacob Faibussowitsch { 422bf7f4e0aSPeter Brune PetscFunctionBegin; 423bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 4246b2b7091SBarry Smith if (linesearch->ops->precheck) { 425dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx); 42638bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed, 4); 427bf7f4e0aSPeter Brune } 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429bf7f4e0aSPeter Brune } 430bf7f4e0aSPeter Brune 431f40b411bSPeter Brune /*@ 432ef46b1a6SStefano Zampini SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 433f40b411bSPeter Brune 434c3339decSBarry Smith Logically Collective 435f40b411bSPeter Brune 436f40b411bSPeter Brune Input Parameters: 4377b1df9c1SPeter Brune + linesearch - The linesearch context 4387b1df9c1SPeter Brune . X - The last solution 4397b1df9c1SPeter Brune . Y - The step direction 4407b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 441f40b411bSPeter Brune 442f40b411bSPeter Brune Output Parameters: 44378bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 44478bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 445f40b411bSPeter Brune 44620f4b53cSBarry Smith Level: developer 44720f4b53cSBarry Smith 448f6dfbefdSBarry Smith Note: 449f6dfbefdSBarry Smith This calls any function provided with `SNESLineSearchSetPreCheck()` 450f6dfbefdSBarry Smith 451db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 452f40b411bSPeter Brune @*/ 453d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 454d71ae5a4SJacob Faibussowitsch { 455bf7f4e0aSPeter Brune PetscFunctionBegin; 456bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 457bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4586b2b7091SBarry Smith if (linesearch->ops->postcheck) { 459dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx); 46038bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5); 46138bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6); 46286d74e61SPeter Brune } 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46486d74e61SPeter Brune } 46586d74e61SPeter Brune 46686d74e61SPeter Brune /*@C 46786d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 46886d74e61SPeter Brune 469c3339decSBarry Smith Logically Collective 47086d74e61SPeter Brune 4714165533cSJose E. Roman Input Parameters: 47286d74e61SPeter Brune + linesearch - linesearch context 47386d74e61SPeter Brune . X - base state for this step 474907376e6SBarry Smith - ctx - context for this function 47586d74e61SPeter Brune 47697bb3fdcSJose E. Roman Input/Output Parameter: 47797bb3fdcSJose E. Roman . Y - correction, possibly modified 47897bb3fdcSJose E. Roman 47997bb3fdcSJose E. Roman Output Parameter: 48097bb3fdcSJose E. Roman . changed - flag indicating that Y was modified 48186d74e61SPeter Brune 48286d74e61SPeter Brune Options Database Key: 483cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 484cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 48586d74e61SPeter Brune 48686d74e61SPeter Brune Level: advanced 48786d74e61SPeter Brune 48886d74e61SPeter Brune Notes: 489f6dfbefdSBarry Smith This function should be passed to `SNESLineSearchSetPreCheck()` 49086d74e61SPeter Brune 49186d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 49286d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 49386d74e61SPeter Brune is generally not useful when using a Newton linearization. 49486d74e61SPeter Brune 495e4094ef1SJacob Faibussowitsch References: 496f6dfbefdSBarry Smith . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 49786d74e61SPeter Brune 498f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()` 49986d74e61SPeter Brune @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) 501d71ae5a4SJacob Faibussowitsch { 50286d74e61SPeter Brune PetscReal angle = *(PetscReal *)linesearch->precheckctx; 50386d74e61SPeter Brune Vec Ylast; 50486d74e61SPeter Brune PetscScalar dot; 50586d74e61SPeter Brune PetscInt iter; 50686d74e61SPeter Brune PetscReal ynorm, ylastnorm, theta, angle_radians; 50786d74e61SPeter Brune SNES snes; 50886d74e61SPeter Brune 50986d74e61SPeter Brune PetscFunctionBegin; 5109566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 5119566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast)); 51286d74e61SPeter Brune if (!Ylast) { 5139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Ylast)); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Ylast)); 51686d74e61SPeter Brune } 5179566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 51886d74e61SPeter Brune if (iter < 2) { 5199566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 52086d74e61SPeter Brune *changed = PETSC_FALSE; 5213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52286d74e61SPeter Brune } 52386d74e61SPeter Brune 5249566063dSJacob Faibussowitsch PetscCall(VecDot(Y, Ylast, &dot)); 5259566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); 5269566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 52786d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 528255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0)); 52986d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 53086d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 53186d74e61SPeter Brune /* Modify the step Y */ 53286d74e61SPeter Brune PetscReal alpha, ydiffnorm; 5339566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ylast, -1.0, Y)); 5349566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm)); 535f85e2ce2SBarry Smith alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0; 5369566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 5379566063dSJacob Faibussowitsch PetscCall(VecScale(Y, alpha)); 5389566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n", (double)(theta * 180. / PETSC_PI), (double)angle, (double)alpha)); 539a47ec85fSBarry Smith *changed = PETSC_TRUE; 54086d74e61SPeter Brune } else { 5419566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle)); 5429566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, Ylast)); 54386d74e61SPeter Brune *changed = PETSC_FALSE; 544bf7f4e0aSPeter Brune } 5453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 546bf7f4e0aSPeter Brune } 547bf7f4e0aSPeter Brune 548f40b411bSPeter Brune /*@ 549cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 550f40b411bSPeter Brune 551c3339decSBarry Smith Collective 552f40b411bSPeter Brune 553f40b411bSPeter Brune Input Parameters: 5548e557f58SPeter Brune + linesearch - The linesearch context 5558e557f58SPeter Brune - Y - The search direction 556f40b411bSPeter Brune 5576b867d5aSJose E. Roman Input/Output Parameters: 5586b867d5aSJose E. Roman + X - The current solution, on output the new solution 5596b867d5aSJose E. Roman . F - The current function, on output the new function 5606b867d5aSJose E. Roman - fnorm - The current norm, on output the new function norm 561f40b411bSPeter Brune 562cd7522eaSPeter Brune Options Database Keys: 5630b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell 564dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 5651fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 5661fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 5673c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 5683c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 569cd7522eaSPeter Brune 570e4094ef1SJacob Faibussowitsch Level: intermediate 57120f4b53cSBarry Smith 572cd7522eaSPeter Brune Notes: 573f6dfbefdSBarry Smith This is typically called from within a `SNESSolve()` implementation in order to 574f6dfbefdSBarry Smith help with convergence of the nonlinear method. Various `SNES` types use line searches 575cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 576cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 577f6dfbefdSBarry Smith application of the line search may invoke `SNESComputeFunction()` several times, and 578cd7522eaSPeter Brune therefore may be fairly expensive. 579cd7522eaSPeter Brune 580f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`, 581db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetType()` 582f40b411bSPeter Brune @*/ 583d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) 584d71ae5a4SJacob Faibussowitsch { 585bf388a1fSBarry Smith PetscFunctionBegin; 586f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 587bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 588bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 589064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Y, VEC_CLASSID, 5); 590bf7f4e0aSPeter Brune 591422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 592bf7f4e0aSPeter Brune 593bf7f4e0aSPeter Brune linesearch->vec_sol = X; 594bf7f4e0aSPeter Brune linesearch->vec_update = Y; 595bf7f4e0aSPeter Brune linesearch->vec_func = F; 596bf7f4e0aSPeter Brune 5979566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(linesearch)); 598bf7f4e0aSPeter Brune 599f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 600bf7f4e0aSPeter Brune 601f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 60248a46eb9SPierre Jolivet else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 603bf7f4e0aSPeter Brune 6049566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 605bf7f4e0aSPeter Brune 606dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, apply); 607bf7f4e0aSPeter Brune 6089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 609bf7f4e0aSPeter Brune 610f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 612bf7f4e0aSPeter Brune } 613bf7f4e0aSPeter Brune 614f40b411bSPeter Brune /*@ 615f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 616f40b411bSPeter Brune 617c3339decSBarry Smith Collective 618f40b411bSPeter Brune 6192fe279fdSBarry Smith Input Parameter: 6208e557f58SPeter Brune . linesearch - The linesearch context 621f40b411bSPeter Brune 62284238204SBarry Smith Level: developer 623f40b411bSPeter Brune 624f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()` 625f40b411bSPeter Brune @*/ 626d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) 627d71ae5a4SJacob Faibussowitsch { 628bf7f4e0aSPeter Brune PetscFunctionBegin; 6293ba16761SJacob Faibussowitsch if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS); 630f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1); 6319371c9d4SSatish Balay if (--((PetscObject)(*linesearch))->refct > 0) { 6329371c9d4SSatish Balay *linesearch = NULL; 6333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6349371c9d4SSatish Balay } 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch)); 6369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset(*linesearch)); 6373ba16761SJacob Faibussowitsch PetscTryTypeMethod(*linesearch, destroy); 6389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 6399566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorCancel((*linesearch))); 6409566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(linesearch)); 6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 642bf7f4e0aSPeter Brune } 643bf7f4e0aSPeter Brune 644f40b411bSPeter Brune /*@ 645dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 646bf7f4e0aSPeter Brune 647c3339decSBarry Smith Logically Collective 648f6dfbefdSBarry Smith 649bf7f4e0aSPeter Brune Input Parameters: 650dcf2fd19SBarry Smith + linesearch - the linesearch object 65120f4b53cSBarry Smith - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor 652bf7f4e0aSPeter Brune 653f6dfbefdSBarry Smith Options Database Key: 654dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 655bf7f4e0aSPeter Brune 656bf7f4e0aSPeter Brune Level: intermediate 657bf7f4e0aSPeter Brune 658e4094ef1SJacob Faibussowitsch Developer Notes: 659f6dfbefdSBarry Smith This monitor is implemented differently than the other line search monitors that are set with 660f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the 661d12e167eSBarry Smith line search that are not visible to the other monitors. 662dcf2fd19SBarry Smith 663f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`, 664f6dfbefdSBarry Smith `SNESLineSearchMonitorSetFromOptions()` 665bf7f4e0aSPeter Brune @*/ 666d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 667d71ae5a4SJacob Faibussowitsch { 668bf7f4e0aSPeter Brune PetscFunctionBegin; 6699566063dSJacob Faibussowitsch if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer)); 6709566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&linesearch->monitor)); 671dcf2fd19SBarry Smith linesearch->monitor = viewer; 6723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 673bf7f4e0aSPeter Brune } 674bf7f4e0aSPeter Brune 675f40b411bSPeter Brune /*@ 676f6dfbefdSBarry Smith SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor. 677f6dfbefdSBarry Smith 678c3339decSBarry Smith Logically Collective 6796a388c36SPeter Brune 680f190f2fcSBarry Smith Input Parameter: 6818e557f58SPeter Brune . linesearch - linesearch context 682f40b411bSPeter Brune 683f190f2fcSBarry Smith Output Parameter: 6848e557f58SPeter Brune . monitor - monitor context 685f40b411bSPeter Brune 686f40b411bSPeter Brune Level: intermediate 687f40b411bSPeter Brune 688f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 689f40b411bSPeter Brune @*/ 690d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 691d71ae5a4SJacob Faibussowitsch { 6926a388c36SPeter Brune PetscFunctionBegin; 693f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 6946a388c36SPeter Brune *monitor = linesearch->monitor; 6953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6966a388c36SPeter Brune } 6976a388c36SPeter Brune 698dcf2fd19SBarry Smith /*@C 699dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 700dcf2fd19SBarry Smith 701c3339decSBarry Smith Collective 702dcf2fd19SBarry Smith 703dcf2fd19SBarry Smith Input Parameters: 704f6dfbefdSBarry Smith + ls - `SNESLineSearch` object you wish to monitor 705dcf2fd19SBarry Smith . name - the monitor type one is seeking 706dcf2fd19SBarry Smith . help - message indicating what monitoring is done 707dcf2fd19SBarry Smith . manual - manual page for the monitor 708dcf2fd19SBarry Smith . monitor - the monitor function 709f6dfbefdSBarry 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` 710dcf2fd19SBarry Smith 711dcf2fd19SBarry Smith Level: developer 712dcf2fd19SBarry Smith 713f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 714db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 715e4094ef1SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 716db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 717c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 718db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 719db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 720dcf2fd19SBarry Smith @*/ 721d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *)) 722d71ae5a4SJacob Faibussowitsch { 723dcf2fd19SBarry Smith PetscViewer viewer; 724dcf2fd19SBarry Smith PetscViewerFormat format; 725dcf2fd19SBarry Smith PetscBool flg; 726dcf2fd19SBarry Smith 727dcf2fd19SBarry Smith PetscFunctionBegin; 7289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg)); 729dcf2fd19SBarry Smith if (flg) { 730d12e167eSBarry Smith PetscViewerAndFormat *vf; 7319566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 7329566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 7331baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls, vf)); 7349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 735dcf2fd19SBarry Smith } 7363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 737dcf2fd19SBarry Smith } 738dcf2fd19SBarry Smith 739f40b411bSPeter Brune /*@ 740f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 741f40b411bSPeter Brune 742c3339decSBarry Smith Logically Collective 743f6dfbefdSBarry Smith 744f899ff85SJose E. Roman Input Parameter: 7458e557f58SPeter Brune . linesearch - linesearch context 746f40b411bSPeter Brune 747f40b411bSPeter Brune Options Database Keys: 7480b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 7493c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 750f6dfbefdSBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`) 75171eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 7521a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 7531a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 7541a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 7551a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 7561a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 757dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 758dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 7598e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 760cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 7611a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 762d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 763f40b411bSPeter Brune 764f40b411bSPeter Brune Level: intermediate 765f40b411bSPeter Brune 766f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 767db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 768f40b411bSPeter Brune @*/ 769d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 770d71ae5a4SJacob Faibussowitsch { 7711a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 772bf7f4e0aSPeter Brune char type[256]; 773bf7f4e0aSPeter Brune PetscBool flg, set; 774dcf2fd19SBarry Smith PetscViewer viewer; 775bf388a1fSBarry Smith 776bf7f4e0aSPeter Brune PetscFunctionBegin; 7779566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll()); 778bf7f4e0aSPeter Brune 779d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch); 780f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 7819566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg)); 782bf7f4e0aSPeter Brune if (flg) { 7839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, type)); 784bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 7859566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, deft)); 786bf7f4e0aSPeter Brune } 787bf7f4e0aSPeter Brune 7889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set)); 789dcf2fd19SBarry Smith if (set) { 7909566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer)); 7919566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 792dcf2fd19SBarry Smith } 7939566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL)); 794bf7f4e0aSPeter Brune 7951a9b3a06SPeter Brune /* tolerances */ 7969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL)); 7979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL)); 7989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL)); 7999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL)); 8009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL)); 8019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL)); 8027a35526eSPeter Brune 8031a9b3a06SPeter Brune /* damping parameters */ 8049566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 8051a9b3a06SPeter Brune 8069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 8071a9b3a06SPeter Brune 8081a9b3a06SPeter Brune /* precheck */ 8099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 8107a35526eSPeter Brune if (set) { 8117a35526eSPeter Brune if (flg) { 8127a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 813f5af7f23SKarl Rupp 814d0609cedSBarry 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)); 8159566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 8167a35526eSPeter Brune } else { 8179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 8187a35526eSPeter Brune } 8197a35526eSPeter Brune } 8209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 8219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 8227a35526eSPeter Brune 823dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 8245a9b6599SPeter Brune 825dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 826d0609cedSBarry Smith PetscOptionsEnd(); 8273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 828bf7f4e0aSPeter Brune } 829bf7f4e0aSPeter Brune 830f40b411bSPeter Brune /*@ 831f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 832f40b411bSPeter Brune 83320f4b53cSBarry Smith Logically Collective 83420f4b53cSBarry Smith 835f40b411bSPeter Brune Input Parameters: 8362fe279fdSBarry Smith + linesearch - linesearch context 8372fe279fdSBarry Smith - viewer - the viewer to display the line search information 838f40b411bSPeter Brune 839f40b411bSPeter Brune Level: intermediate 840f40b411bSPeter Brune 841f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 842f40b411bSPeter Brune @*/ 843d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 844d71ae5a4SJacob Faibussowitsch { 8457f1410a3SPeter Brune PetscBool iascii; 846bf388a1fSBarry Smith 847bf7f4e0aSPeter Brune PetscFunctionBegin; 8487f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 84948a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 8507f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 8517f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2); 852f40b411bSPeter Brune 8539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 8547f1410a3SPeter Brune if (iascii) { 8559566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 8569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 857dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer); 8589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 8599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 8609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 86163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 8626b2b7091SBarry Smith if (linesearch->ops->precheck) { 8636b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 86463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 8657f1410a3SPeter Brune } else { 86663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 8677f1410a3SPeter Brune } 8687f1410a3SPeter Brune } 86948a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 8707f1410a3SPeter Brune } 8713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 872bf7f4e0aSPeter Brune } 873bf7f4e0aSPeter Brune 874ea5d4fccSPeter Brune /*@C 875a80ff896SJed Brown SNESLineSearchGetType - Gets the linesearch type 876a80ff896SJed Brown 877c3339decSBarry Smith Logically Collective 878a80ff896SJed Brown 8792fe279fdSBarry Smith Input Parameter: 880a80ff896SJed Brown . linesearch - linesearch context 881a80ff896SJed Brown 8822fe279fdSBarry Smith Output Parameter: 8832fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set 884a80ff896SJed Brown 885a80ff896SJed Brown Level: intermediate 886a80ff896SJed Brown 887e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 888a80ff896SJed Brown @*/ 889d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 890d71ae5a4SJacob Faibussowitsch { 891a80ff896SJed Brown PetscFunctionBegin; 892a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 8934f572ea9SToby Isaac PetscAssertPointer(type, 2); 894a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 8953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 896a80ff896SJed Brown } 897a80ff896SJed Brown 898a80ff896SJed Brown /*@C 899f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 900f40b411bSPeter Brune 901c3339decSBarry Smith Logically Collective 902f190f2fcSBarry Smith 903f40b411bSPeter Brune Input Parameters: 9048e557f58SPeter Brune + linesearch - linesearch context 905ceaaa498SBarry Smith - type - The type of line search to be used, see `SNESLineSearchType` 9061fe24845SBarry Smith 9073c7db156SBarry Smith Options Database Key: 9080b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 909cd7522eaSPeter Brune 910f40b411bSPeter Brune Level: intermediate 911f40b411bSPeter Brune 912e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 913f40b411bSPeter Brune @*/ 914d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 915d71ae5a4SJacob Faibussowitsch { 916bf7f4e0aSPeter Brune PetscBool match; 9175f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 918bf7f4e0aSPeter Brune 919bf7f4e0aSPeter Brune PetscFunctionBegin; 920f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9214f572ea9SToby Isaac PetscAssertPointer(type, 2); 922bf7f4e0aSPeter Brune 9239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 9243ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 925bf7f4e0aSPeter Brune 9269566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 9276adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 928bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 929bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 9309566063dSJacob Faibussowitsch PetscCall((*(linesearch)->ops->destroy)(linesearch)); 9310298fd71SBarry Smith linesearch->ops->destroy = NULL; 932bf7f4e0aSPeter Brune } 933f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9349e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9359e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9369e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9379e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 938bf7f4e0aSPeter Brune 9399566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 9409566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 9413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 942bf7f4e0aSPeter Brune } 943bf7f4e0aSPeter Brune 944f40b411bSPeter Brune /*@ 945f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 946f40b411bSPeter Brune 947f40b411bSPeter Brune Input Parameters: 9488e557f58SPeter Brune + linesearch - linesearch context 949f40b411bSPeter Brune - snes - The snes instance 950f40b411bSPeter Brune 95178bcb3b5SPeter Brune Level: developer 95278bcb3b5SPeter Brune 953f6dfbefdSBarry Smith Note: 954f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 955f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 95678bcb3b5SPeter Brune implementations. 957f40b411bSPeter Brune 958f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 959f40b411bSPeter Brune @*/ 960d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 961d71ae5a4SJacob Faibussowitsch { 962bf7f4e0aSPeter Brune PetscFunctionBegin; 963f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 964bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 965bf7f4e0aSPeter Brune linesearch->snes = snes; 9663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 967bf7f4e0aSPeter Brune } 968bf7f4e0aSPeter Brune 969f40b411bSPeter Brune /*@ 970f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 971f6dfbefdSBarry Smith 972f6dfbefdSBarry Smith Not Collective 973f40b411bSPeter Brune 9742fe279fdSBarry Smith Input Parameter: 9758e557f58SPeter Brune . linesearch - linesearch context 976f40b411bSPeter Brune 9772fe279fdSBarry Smith Output Parameter: 9782fe279fdSBarry Smith . snes - The `SNES` instance 979f40b411bSPeter Brune 9808141a3b9SPeter Brune Level: developer 981f40b411bSPeter Brune 98242747ad1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`, `SNES` 983f40b411bSPeter Brune @*/ 984d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 985d71ae5a4SJacob Faibussowitsch { 986bf7f4e0aSPeter Brune PetscFunctionBegin; 987f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 9884f572ea9SToby Isaac PetscAssertPointer(snes, 2); 989bf7f4e0aSPeter Brune *snes = linesearch->snes; 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 991bf7f4e0aSPeter Brune } 992bf7f4e0aSPeter Brune 993f40b411bSPeter Brune /*@ 994f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 995f40b411bSPeter Brune 996f6dfbefdSBarry Smith Not Collective 997f6dfbefdSBarry Smith 9982fe279fdSBarry Smith Input Parameter: 9998e557f58SPeter Brune . linesearch - linesearch context 1000f40b411bSPeter Brune 10012fe279fdSBarry Smith Output Parameter: 1002f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()` 1003f40b411bSPeter Brune 100478bcb3b5SPeter Brune Level: advanced 100578bcb3b5SPeter Brune 1006f6dfbefdSBarry Smith Note: 10078e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 100878bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 1009f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the 101078bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 101178bcb3b5SPeter Brune 1012f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1013f40b411bSPeter Brune @*/ 1014d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1015d71ae5a4SJacob Faibussowitsch { 10166a388c36SPeter Brune PetscFunctionBegin; 1017f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10184f572ea9SToby Isaac PetscAssertPointer(lambda, 2); 10196a388c36SPeter Brune *lambda = linesearch->lambda; 10203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10216a388c36SPeter Brune } 10226a388c36SPeter Brune 1023f40b411bSPeter Brune /*@ 1024f6dfbefdSBarry Smith SNESLineSearchSetLambda - Sets the linesearch steplength 1025f40b411bSPeter Brune 1026f40b411bSPeter Brune Input Parameters: 10278e557f58SPeter Brune + linesearch - linesearch context 1028f40b411bSPeter Brune - lambda - The last steplength. 1029f40b411bSPeter Brune 103020f4b53cSBarry Smith Level: advanced 103120f4b53cSBarry Smith 1032f6dfbefdSBarry Smith Note: 1033f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()` 1034f6dfbefdSBarry Smith to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1035cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1036cd7522eaSPeter Brune as an inner scaling parameter. 1037cd7522eaSPeter Brune 1038f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()` 1039f40b411bSPeter Brune @*/ 1040d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1041d71ae5a4SJacob Faibussowitsch { 10426a388c36SPeter Brune PetscFunctionBegin; 1043f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 10446a388c36SPeter Brune linesearch->lambda = lambda; 10453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10466a388c36SPeter Brune } 10476a388c36SPeter Brune 1048f40b411bSPeter Brune /*@ 1049ceaaa498SBarry Smith SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. 1050f40b411bSPeter Brune 1051f6dfbefdSBarry Smith Not Collective 1052f6dfbefdSBarry Smith 1053f899ff85SJose E. Roman Input Parameter: 10548e557f58SPeter Brune . linesearch - linesearch context 1055f40b411bSPeter Brune 1056f40b411bSPeter Brune Output Parameters: 1057b13b64c2SBarry Smith + steptol - The minimum steplength 10586cc8e53bSPeter Brune . maxstep - The maximum steplength 1059516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1060516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1061516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1062e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search 1063f40b411bSPeter Brune 106478bcb3b5SPeter Brune Level: intermediate 106578bcb3b5SPeter Brune 1066f6dfbefdSBarry Smith Note: 106778bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 10683c7d6663SPeter Brune the type requires. 1069516fe3c3SPeter Brune 1070f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1071f40b411bSPeter Brune @*/ 1072d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1073d71ae5a4SJacob Faibussowitsch { 10746a388c36SPeter Brune PetscFunctionBegin; 1075f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1076516fe3c3SPeter Brune if (steptol) { 10774f572ea9SToby Isaac PetscAssertPointer(steptol, 2); 10786a388c36SPeter Brune *steptol = linesearch->steptol; 1079516fe3c3SPeter Brune } 1080516fe3c3SPeter Brune if (maxstep) { 10814f572ea9SToby Isaac PetscAssertPointer(maxstep, 3); 1082b13b64c2SBarry Smith *maxstep = linesearch->maxstep; 1083516fe3c3SPeter Brune } 1084516fe3c3SPeter Brune if (rtol) { 10854f572ea9SToby Isaac PetscAssertPointer(rtol, 4); 1086516fe3c3SPeter Brune *rtol = linesearch->rtol; 1087516fe3c3SPeter Brune } 1088516fe3c3SPeter Brune if (atol) { 10894f572ea9SToby Isaac PetscAssertPointer(atol, 5); 1090516fe3c3SPeter Brune *atol = linesearch->atol; 1091516fe3c3SPeter Brune } 1092516fe3c3SPeter Brune if (ltol) { 10934f572ea9SToby Isaac PetscAssertPointer(ltol, 6); 1094516fe3c3SPeter Brune *ltol = linesearch->ltol; 1095516fe3c3SPeter Brune } 1096516fe3c3SPeter Brune if (max_its) { 10974f572ea9SToby Isaac PetscAssertPointer(max_its, 7); 1098516fe3c3SPeter Brune *max_its = linesearch->max_its; 1099516fe3c3SPeter Brune } 11003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11016a388c36SPeter Brune } 11026a388c36SPeter Brune 1103f40b411bSPeter Brune /*@ 1104ceaaa498SBarry Smith SNESLineSearchSetTolerances - Sets the tolerances for the linesearch. 1105f40b411bSPeter Brune 1106f6dfbefdSBarry Smith Collective 1107f6dfbefdSBarry Smith 1108f40b411bSPeter Brune Input Parameters: 11098e557f58SPeter Brune + linesearch - linesearch context 1110516fe3c3SPeter Brune . steptol - The minimum steplength 11116cc8e53bSPeter Brune . maxstep - The maximum steplength 1112516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1113516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1114516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1115e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search 1116f40b411bSPeter Brune 111720f4b53cSBarry Smith Level: intermediate 111820f4b53cSBarry Smith 1119f6dfbefdSBarry Smith Note: 1120f6dfbefdSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in 112178bcb3b5SPeter Brune place of an argument. 1122f40b411bSPeter Brune 1123f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1124f40b411bSPeter Brune @*/ 1125d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 1126d71ae5a4SJacob Faibussowitsch { 11276a388c36SPeter Brune PetscFunctionBegin; 1128f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1129d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1130d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1131d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1132d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1133d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1134d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch, max_its, 7); 1135d3952184SSatish Balay 113613bcc0bdSJacob Faibussowitsch if (steptol != (PetscReal)PETSC_DEFAULT) { 11375f80ce2aSJacob Faibussowitsch PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 11386a388c36SPeter Brune linesearch->steptol = steptol; 1139d3952184SSatish Balay } 1140d3952184SSatish Balay 114113bcc0bdSJacob Faibussowitsch if (maxstep != (PetscReal)PETSC_DEFAULT) { 11425f80ce2aSJacob Faibussowitsch PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1143516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1144d3952184SSatish Balay } 1145d3952184SSatish Balay 114613bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) { 11472061ca32SJunchao 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); 1148516fe3c3SPeter Brune linesearch->rtol = rtol; 1149d3952184SSatish Balay } 1150d3952184SSatish Balay 115113bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) { 11525f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1153516fe3c3SPeter Brune linesearch->atol = atol; 1154d3952184SSatish Balay } 1155d3952184SSatish Balay 115613bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) { 11575f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1158516fe3c3SPeter Brune linesearch->ltol = ltol; 1159d3952184SSatish Balay } 1160d3952184SSatish Balay 1161d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 116263a3b9bcSJacob Faibussowitsch PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its); 1163516fe3c3SPeter Brune linesearch->max_its = max_its; 1164d3952184SSatish Balay } 11653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11666a388c36SPeter Brune } 11676a388c36SPeter Brune 1168f40b411bSPeter Brune /*@ 1169f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1170f40b411bSPeter Brune 11712fe279fdSBarry Smith Input Parameter: 11728e557f58SPeter Brune . linesearch - linesearch context 1173f40b411bSPeter Brune 11742fe279fdSBarry Smith Output Parameter: 11758e557f58SPeter Brune . damping - The damping parameter 1176f40b411bSPeter Brune 117778bcb3b5SPeter Brune Level: advanced 1178f40b411bSPeter Brune 1179db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN` 1180f40b411bSPeter Brune @*/ 1181d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1182d71ae5a4SJacob Faibussowitsch { 11836a388c36SPeter Brune PetscFunctionBegin; 1184f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 11854f572ea9SToby Isaac PetscAssertPointer(damping, 2); 11866a388c36SPeter Brune *damping = linesearch->damping; 11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11886a388c36SPeter Brune } 11896a388c36SPeter Brune 1190f40b411bSPeter Brune /*@ 1191fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1192f40b411bSPeter Brune 1193f40b411bSPeter Brune Input Parameters: 119403fd4120SBarry Smith + linesearch - linesearch context 119503fd4120SBarry Smith - damping - The damping parameter 1196f40b411bSPeter Brune 1197f6dfbefdSBarry Smith Options Database Key: 1198f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value 1199f6dfbefdSBarry Smith 1200f40b411bSPeter Brune Level: intermediate 1201f40b411bSPeter Brune 1202f6dfbefdSBarry Smith Note: 1203f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1204cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 120578bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1206cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1207cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1208cd7522eaSPeter Brune 1209f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()` 1210f40b411bSPeter Brune @*/ 1211d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1212d71ae5a4SJacob Faibussowitsch { 12136a388c36SPeter Brune PetscFunctionBegin; 1214f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12156a388c36SPeter Brune linesearch->damping = damping; 12163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12176a388c36SPeter Brune } 12186a388c36SPeter Brune 121959405d5eSPeter Brune /*@ 122059405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 122159405d5eSPeter Brune 1222f6dfbefdSBarry Smith Input Parameter: 122378bcb3b5SPeter Brune . linesearch - linesearch context 122459405d5eSPeter Brune 1225f6dfbefdSBarry Smith Output Parameter: 12268e557f58SPeter Brune . order - The order 122759405d5eSPeter Brune 122859405d5eSPeter Brune Level: intermediate 122959405d5eSPeter Brune 1230f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()` 123159405d5eSPeter Brune @*/ 1232d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1233d71ae5a4SJacob Faibussowitsch { 123459405d5eSPeter Brune PetscFunctionBegin; 123559405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 12364f572ea9SToby Isaac PetscAssertPointer(order, 2); 123759405d5eSPeter Brune *order = linesearch->order; 12383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123959405d5eSPeter Brune } 124059405d5eSPeter Brune 124159405d5eSPeter Brune /*@ 12421f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 124359405d5eSPeter Brune 124459405d5eSPeter Brune Input Parameters: 1245a2b725a8SWilliam Gropp + linesearch - linesearch context 1246ceaaa498SBarry Smith - order - The order 124759405d5eSPeter Brune 124859405d5eSPeter Brune Level: intermediate 124959405d5eSPeter Brune 1250ceaaa498SBarry Smith Values for `order`\: 1251f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1252f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1253f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 125478bcb3b5SPeter Brune 1255ceaaa498SBarry Smith Note: 1256ceaaa498SBarry Smith These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP` 125759405d5eSPeter Brune 1258f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 125959405d5eSPeter Brune @*/ 1260d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1261d71ae5a4SJacob Faibussowitsch { 126259405d5eSPeter Brune PetscFunctionBegin; 126359405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 126459405d5eSPeter Brune linesearch->order = order; 12653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126659405d5eSPeter Brune } 126759405d5eSPeter Brune 1268f40b411bSPeter Brune /*@ 1269f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1270f40b411bSPeter Brune 1271f6dfbefdSBarry Smith Not Collective 1272f6dfbefdSBarry Smith 1273f899ff85SJose E. Roman Input Parameter: 127478bcb3b5SPeter Brune . linesearch - linesearch context 1275f40b411bSPeter Brune 1276f40b411bSPeter Brune Output Parameters: 1277f40b411bSPeter Brune + xnorm - The norm of the current solution 1278a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1279f40b411bSPeter Brune - ynorm - The norm of the current update 1280f40b411bSPeter Brune 128178bcb3b5SPeter Brune Level: developer 1282f40b411bSPeter Brune 1283a68bbae5SBarry Smith Note: 1284a68bbae5SBarry Smith Some values may not be up to date at particular points in the code. 1285a68bbae5SBarry Smith 1286a68bbae5SBarry Smith This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1287a68bbae5SBarry Smith computed values. 1288a68bbae5SBarry Smith 1289f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1290f40b411bSPeter Brune @*/ 1291d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1292d71ae5a4SJacob Faibussowitsch { 1293bf7f4e0aSPeter Brune PetscFunctionBegin; 1294f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1295f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1296f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1297f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 12983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1299bf7f4e0aSPeter Brune } 1300bf7f4e0aSPeter Brune 1301f40b411bSPeter Brune /*@ 1302f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1303f40b411bSPeter Brune 1304c3339decSBarry Smith Collective 1305f6dfbefdSBarry Smith 1306f40b411bSPeter Brune Input Parameters: 130778bcb3b5SPeter Brune + linesearch - linesearch context 1308f40b411bSPeter Brune . xnorm - The norm of the current solution 1309a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1310f40b411bSPeter Brune - ynorm - The norm of the current update 1311f40b411bSPeter Brune 1312f6dfbefdSBarry Smith Level: developer 1313f40b411bSPeter Brune 1314f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1315f40b411bSPeter Brune @*/ 1316d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1317d71ae5a4SJacob Faibussowitsch { 13186a388c36SPeter Brune PetscFunctionBegin; 1319f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 13206a388c36SPeter Brune linesearch->xnorm = xnorm; 13216a388c36SPeter Brune linesearch->fnorm = fnorm; 13226a388c36SPeter Brune linesearch->ynorm = ynorm; 13233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13246a388c36SPeter Brune } 13256a388c36SPeter Brune 1326f40b411bSPeter Brune /*@ 1327f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1328f40b411bSPeter Brune 13292fe279fdSBarry Smith Input Parameter: 133078bcb3b5SPeter Brune . linesearch - linesearch context 1331f40b411bSPeter Brune 133220f4b53cSBarry Smith Options Database Key: 13338e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1334f40b411bSPeter Brune 1335f40b411bSPeter Brune Level: intermediate 1336f40b411bSPeter Brune 1337f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1338f40b411bSPeter Brune @*/ 1339d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1340d71ae5a4SJacob Faibussowitsch { 13419bd66eb0SPeter Brune SNES snes; 1342bf388a1fSBarry Smith 13436a388c36SPeter Brune PetscFunctionBegin; 13446a388c36SPeter Brune if (linesearch->norms) { 13459bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 13469566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 13479566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13489566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13499566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 13509bd66eb0SPeter Brune } else { 13519566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13529566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13539566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13549566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13559566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13569566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13576a388c36SPeter Brune } 13589bd66eb0SPeter Brune } 13593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13606a388c36SPeter Brune } 13616a388c36SPeter Brune 13626f263ca3SPeter Brune /*@ 13636f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 13646f263ca3SPeter Brune 13656f263ca3SPeter Brune Input Parameters: 136678bcb3b5SPeter Brune + linesearch - linesearch context 136778bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 13686f263ca3SPeter Brune 136920f4b53cSBarry Smith Options Database Key: 13700b00b554SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch 13716f263ca3SPeter Brune 137220f4b53cSBarry Smith Level: intermediate 137320f4b53cSBarry Smith 1374f6dfbefdSBarry Smith Note: 1375f6dfbefdSBarry 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. 13766f263ca3SPeter Brune 1377f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 13786f263ca3SPeter Brune @*/ 1379d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1380d71ae5a4SJacob Faibussowitsch { 13816f263ca3SPeter Brune PetscFunctionBegin; 13826f263ca3SPeter Brune linesearch->norms = flg; 13833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13846f263ca3SPeter Brune } 13856f263ca3SPeter Brune 1386f40b411bSPeter Brune /*@ 1387f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1388f6dfbefdSBarry Smith 1389f6dfbefdSBarry Smith Not Collective on but the vectors are parallel 1390f40b411bSPeter Brune 1391f899ff85SJose E. Roman Input Parameter: 139278bcb3b5SPeter Brune . linesearch - linesearch context 1393f40b411bSPeter Brune 1394f40b411bSPeter Brune Output Parameters: 13956232e825SPeter Brune + X - Solution vector 13966232e825SPeter Brune . F - Function vector 13976232e825SPeter Brune . Y - Search direction vector 13986232e825SPeter Brune . W - Solution work vector 13996232e825SPeter Brune - G - Function work vector 14006232e825SPeter Brune 140120f4b53cSBarry Smith Level: advanced 140220f4b53cSBarry Smith 14037bba9028SPeter Brune Notes: 140420f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a 140520f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the 140620f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the 14076232e825SPeter Brune function evaluated at the new solution. 1408f40b411bSPeter Brune 1409f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 14102a7a6963SBarry Smith 1411f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1412f40b411bSPeter Brune @*/ 1413d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1414d71ae5a4SJacob Faibussowitsch { 14156a388c36SPeter Brune PetscFunctionBegin; 1416f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14176a388c36SPeter Brune if (X) { 14184f572ea9SToby Isaac PetscAssertPointer(X, 2); 14196a388c36SPeter Brune *X = linesearch->vec_sol; 14206a388c36SPeter Brune } 14216a388c36SPeter Brune if (F) { 14224f572ea9SToby Isaac PetscAssertPointer(F, 3); 14236a388c36SPeter Brune *F = linesearch->vec_func; 14246a388c36SPeter Brune } 14256a388c36SPeter Brune if (Y) { 14264f572ea9SToby Isaac PetscAssertPointer(Y, 4); 14276a388c36SPeter Brune *Y = linesearch->vec_update; 14286a388c36SPeter Brune } 14296a388c36SPeter Brune if (W) { 14304f572ea9SToby Isaac PetscAssertPointer(W, 5); 14316a388c36SPeter Brune *W = linesearch->vec_sol_new; 14326a388c36SPeter Brune } 14336a388c36SPeter Brune if (G) { 14344f572ea9SToby Isaac PetscAssertPointer(G, 6); 14356a388c36SPeter Brune *G = linesearch->vec_func_new; 14366a388c36SPeter Brune } 14373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14386a388c36SPeter Brune } 14396a388c36SPeter Brune 1440f40b411bSPeter Brune /*@ 1441f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1442f6dfbefdSBarry Smith 1443c3339decSBarry Smith Logically Collective 1444f40b411bSPeter Brune 1445f40b411bSPeter Brune Input Parameters: 144678bcb3b5SPeter Brune + linesearch - linesearch context 14476232e825SPeter Brune . X - Solution vector 14486232e825SPeter Brune . F - Function vector 14496232e825SPeter Brune . Y - Search direction vector 14506232e825SPeter Brune . W - Solution work vector 14516232e825SPeter Brune - G - Function work vector 1452f40b411bSPeter Brune 145378bcb3b5SPeter Brune Level: advanced 1454f40b411bSPeter Brune 1455f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1456f40b411bSPeter Brune @*/ 1457d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1458d71ae5a4SJacob Faibussowitsch { 14596a388c36SPeter Brune PetscFunctionBegin; 1460f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 14616a388c36SPeter Brune if (X) { 14626a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 14636a388c36SPeter Brune linesearch->vec_sol = X; 14646a388c36SPeter Brune } 14656a388c36SPeter Brune if (F) { 14666a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 14676a388c36SPeter Brune linesearch->vec_func = F; 14686a388c36SPeter Brune } 14696a388c36SPeter Brune if (Y) { 14706a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 14716a388c36SPeter Brune linesearch->vec_update = Y; 14726a388c36SPeter Brune } 14736a388c36SPeter Brune if (W) { 14746a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 14756a388c36SPeter Brune linesearch->vec_sol_new = W; 14766a388c36SPeter Brune } 14776a388c36SPeter Brune if (G) { 14786a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 14796a388c36SPeter Brune linesearch->vec_func_new = G; 14806a388c36SPeter Brune } 14813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14826a388c36SPeter Brune } 14836a388c36SPeter Brune 1484e7058c64SPeter Brune /*@C 1485f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1486f6dfbefdSBarry Smith `SNESLineSearch` options in the database. 1487e7058c64SPeter Brune 1488c3339decSBarry Smith Logically Collective 1489e7058c64SPeter Brune 1490e7058c64SPeter Brune Input Parameters: 1491f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1492e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1493e7058c64SPeter Brune 149420f4b53cSBarry Smith Level: advanced 149520f4b53cSBarry Smith 1496f6dfbefdSBarry Smith Note: 1497e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1498e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1499e7058c64SPeter Brune 1500f6dfbefdSBarry Smith .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1501e7058c64SPeter Brune @*/ 1502d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1503d71ae5a4SJacob Faibussowitsch { 1504e7058c64SPeter Brune PetscFunctionBegin; 1505f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15069566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 15073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1508e7058c64SPeter Brune } 1509e7058c64SPeter Brune 1510e7058c64SPeter Brune /*@C 1511f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1512f1c6b773SPeter Brune SNESLineSearch options in the database. 1513e7058c64SPeter Brune 1514e7058c64SPeter Brune Not Collective 1515e7058c64SPeter Brune 1516e7058c64SPeter Brune Input Parameter: 1517f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context 1518e7058c64SPeter Brune 1519e7058c64SPeter Brune Output Parameter: 1520e7058c64SPeter Brune . prefix - pointer to the prefix string used 1521e7058c64SPeter Brune 1522e7058c64SPeter Brune Level: advanced 1523e7058c64SPeter Brune 1524e4094ef1SJacob Faibussowitsch Fortran Notes: 152520f4b53cSBarry Smith The user should pass in a string 'prefix' of 152620f4b53cSBarry Smith sufficient length to hold the prefix. 152720f4b53cSBarry Smith 1528f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1529e7058c64SPeter Brune @*/ 1530d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1531d71ae5a4SJacob Faibussowitsch { 1532e7058c64SPeter Brune PetscFunctionBegin; 1533f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 15353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1536e7058c64SPeter Brune } 1537bf7f4e0aSPeter Brune 15388d359177SBarry Smith /*@C 1539f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1540f40b411bSPeter Brune 1541d8d19677SJose E. Roman Input Parameters: 1542f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context 1543f40b411bSPeter Brune - nwork - the number of work vectors 1544f40b411bSPeter Brune 1545f40b411bSPeter Brune Level: developer 1546f40b411bSPeter Brune 1547f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetWorkVecs()` 1548f40b411bSPeter Brune @*/ 1549d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1550d71ae5a4SJacob Faibussowitsch { 1551bf7f4e0aSPeter Brune PetscFunctionBegin; 15520fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 15539566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 15543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1555bf7f4e0aSPeter Brune } 1556bf7f4e0aSPeter Brune 1557f40b411bSPeter Brune /*@ 1558422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1559f40b411bSPeter Brune 15602fe279fdSBarry Smith Input Parameter: 156178bcb3b5SPeter Brune . linesearch - linesearch context 1562f40b411bSPeter Brune 15632fe279fdSBarry Smith Output Parameter: 1564422a814eSBarry Smith . result - The success or failure status 1565f40b411bSPeter Brune 156620f4b53cSBarry Smith Level: developer 156720f4b53cSBarry Smith 1568f6dfbefdSBarry Smith Note: 1569f6dfbefdSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed 1570cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1571cd7522eaSPeter Brune 1572f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1573f40b411bSPeter Brune @*/ 1574d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1575d71ae5a4SJacob Faibussowitsch { 1576bf7f4e0aSPeter Brune PetscFunctionBegin; 1577f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 15784f572ea9SToby Isaac PetscAssertPointer(result, 2); 1579422a814eSBarry Smith *result = linesearch->result; 15803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1581bf7f4e0aSPeter Brune } 1582bf7f4e0aSPeter Brune 1583f40b411bSPeter Brune /*@ 1584422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1585f40b411bSPeter Brune 1586f40b411bSPeter Brune Input Parameters: 158778bcb3b5SPeter Brune + linesearch - linesearch context 1588422a814eSBarry Smith - result - The success or failure status 1589f40b411bSPeter Brune 159020f4b53cSBarry Smith Level: developer 159120f4b53cSBarry Smith 1592f6dfbefdSBarry Smith Note: 1593f6dfbefdSBarry Smith This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1594cd7522eaSPeter Brune the success or failure of the line search method. 1595cd7522eaSPeter Brune 1596f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1597f40b411bSPeter Brune @*/ 1598d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1599d71ae5a4SJacob Faibussowitsch { 16006a388c36SPeter Brune PetscFunctionBegin; 16015fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1602422a814eSBarry Smith linesearch->result = result; 16033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16046a388c36SPeter Brune } 16056a388c36SPeter Brune 1606ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 16079bd66eb0SPeter Brune /*@C 1608f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16099bd66eb0SPeter Brune 1610c3339decSBarry Smith Logically Collective 1611f6dfbefdSBarry Smith 16129bd66eb0SPeter Brune Input Parameters: 1613ceaaa498SBarry Smith + linesearch - the linesearch object 16149bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 16159bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16169bd66eb0SPeter Brune 161720f4b53cSBarry Smith Calling sequence of `projectfunc`: 1618ceaaa498SBarry Smith $ PetscErrorCode projectfunc(SNES snes, Vec X) 16199bd66eb0SPeter Brune + snes - nonlinear context 162020f4b53cSBarry Smith - X - current solution, store the projected solution here 16219bd66eb0SPeter Brune 162220f4b53cSBarry Smith Calling sequence of `normfunc`: 1623ceaaa498SBarry Smith $ PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm) 16249bd66eb0SPeter Brune + snes - nonlinear context 16259bd66eb0SPeter Brune . X - current solution 162620f4b53cSBarry Smith . F - current residual 1627ceaaa498SBarry Smith - fnorm - the computed norm 16289bd66eb0SPeter Brune 1629f6dfbefdSBarry Smith Level: advanced 16309bd66eb0SPeter Brune 1631ceaaa498SBarry Smith Notes: 163220f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 163320f4b53cSBarry Smith 163420f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated 163520f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`. 163620f4b53cSBarry Smith 1637f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 16389bd66eb0SPeter Brune @*/ 1639d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1640d71ae5a4SJacob Faibussowitsch { 16419bd66eb0SPeter Brune PetscFunctionBegin; 1642f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 16439bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 16449bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 16453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16469bd66eb0SPeter Brune } 16479bd66eb0SPeter Brune 16489bd66eb0SPeter Brune /*@C 1649f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 16509bd66eb0SPeter Brune 1651f6dfbefdSBarry Smith Not Collective 1652f6dfbefdSBarry Smith 1653f899ff85SJose E. Roman Input Parameter: 1654f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()` 16559bd66eb0SPeter Brune 16569bd66eb0SPeter Brune Output Parameters: 16579bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16589bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16599bd66eb0SPeter Brune 1660f6dfbefdSBarry Smith Level: advanced 16619bd66eb0SPeter Brune 1662f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 16639bd66eb0SPeter Brune @*/ 1664d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1665d71ae5a4SJacob Faibussowitsch { 16669bd66eb0SPeter Brune PetscFunctionBegin; 16679bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16689bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16709bd66eb0SPeter Brune } 16719bd66eb0SPeter Brune 1672bf7f4e0aSPeter Brune /*@C 1673ceaaa498SBarry Smith SNESLineSearchRegister - register a line search type 1674ceaaa498SBarry Smith 1675ceaaa498SBarry Smith Input Parameters: 1676ceaaa498SBarry Smith + sname - name of the line search 1677ceaaa498SBarry Smith - function - the creation function for that type 1678ceaaa498SBarry Smith 1679ceaaa498SBarry Smith Calling sequence of `function`: 1680ceaaa498SBarry Smith . ls - the linesearch context 1681bf7f4e0aSPeter Brune 1682bf7f4e0aSPeter Brune Level: advanced 1683f6dfbefdSBarry Smith 1684f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1685bf7f4e0aSPeter Brune @*/ 1686ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls)) 1687d71ae5a4SJacob Faibussowitsch { 1688bf7f4e0aSPeter Brune PetscFunctionBegin; 16899566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 16909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 16913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1692bf7f4e0aSPeter Brune } 1693