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 /*@ 10dcf2fd19SBarry Smith SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object. 11dcf2fd19SBarry Smith 12dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 13dcf2fd19SBarry Smith 14dcf2fd19SBarry Smith Input Parameters: 15dcf2fd19SBarry 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 19dcf2fd19SBarry Smith into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those 20dcf2fd19SBarry Smith set via the options database 21dcf2fd19SBarry Smith 22dcf2fd19SBarry Smith Notes: 23dcf2fd19SBarry Smith There is no way to clear one specific monitor from a SNESLineSearch object. 24dcf2fd19SBarry Smith 25dcf2fd19SBarry Smith This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel 26dcf2fd19SBarry Smith that one. 27dcf2fd19SBarry Smith 28dcf2fd19SBarry Smith Level: intermediate 29dcf2fd19SBarry Smith 30db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()` 31dcf2fd19SBarry Smith @*/ 32dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 33dcf2fd19SBarry Smith { 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++) { 39dcf2fd19SBarry Smith if (ls->monitordestroy[i]) { 409566063dSJacob Faibussowitsch PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i])); 41dcf2fd19SBarry Smith } 42dcf2fd19SBarry Smith } 43dcf2fd19SBarry Smith ls->numbermonitors = 0; 44dcf2fd19SBarry Smith PetscFunctionReturn(0); 45dcf2fd19SBarry Smith } 46dcf2fd19SBarry Smith 47dcf2fd19SBarry Smith /*@ 48dcf2fd19SBarry Smith SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 49dcf2fd19SBarry Smith 50dcf2fd19SBarry Smith Collective on SNES 51dcf2fd19SBarry Smith 52dcf2fd19SBarry Smith Input Parameters: 53dcf2fd19SBarry Smith . ls - the linesearch object 54dcf2fd19SBarry Smith 55dcf2fd19SBarry Smith Notes: 56dcf2fd19SBarry Smith This routine is called by the SNES implementations. 57dcf2fd19SBarry Smith It does not typically need to be called by the user. 58dcf2fd19SBarry Smith 59dcf2fd19SBarry Smith Level: developer 60dcf2fd19SBarry Smith 61db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()` 62dcf2fd19SBarry Smith @*/ 63dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 64dcf2fd19SBarry Smith { 65dcf2fd19SBarry Smith PetscInt i,n = ls->numbermonitors; 66dcf2fd19SBarry Smith 67dcf2fd19SBarry Smith PetscFunctionBegin; 68dcf2fd19SBarry Smith for (i=0; i<n; i++) { 699566063dSJacob Faibussowitsch PetscCall((*ls->monitorftns[i])(ls,ls->monitorcontext[i])); 70dcf2fd19SBarry Smith } 71dcf2fd19SBarry Smith PetscFunctionReturn(0); 72dcf2fd19SBarry Smith } 73dcf2fd19SBarry Smith 74dcf2fd19SBarry Smith /*@C 75dcf2fd19SBarry Smith SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 76dcf2fd19SBarry Smith iteration of the nonlinear solver to display the iteration's 77dcf2fd19SBarry Smith progress. 78dcf2fd19SBarry Smith 79dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 80dcf2fd19SBarry Smith 81dcf2fd19SBarry Smith Input Parameters: 82dcf2fd19SBarry Smith + ls - the SNESLineSearch context 83dcf2fd19SBarry Smith . f - the monitor function 84dcf2fd19SBarry Smith . mctx - [optional] user-defined context for private data for the 85dcf2fd19SBarry Smith monitor routine (use NULL if no context is desired) 86dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context 87dcf2fd19SBarry Smith (may be NULL) 88dcf2fd19SBarry Smith 89dcf2fd19SBarry Smith Notes: 90dcf2fd19SBarry Smith Several different monitoring routines may be set by calling 91dcf2fd19SBarry Smith SNESLineSearchMonitorSet() multiple times; all will be called in the 92dcf2fd19SBarry Smith order in which they were set. 93dcf2fd19SBarry Smith 9495452b02SPatrick Sanan Fortran Notes: 9595452b02SPatrick Sanan Only a single monitor function can be set for each SNESLineSearch object 96dcf2fd19SBarry Smith 97dcf2fd19SBarry Smith Level: intermediate 98dcf2fd19SBarry Smith 99db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()` 100dcf2fd19SBarry Smith @*/ 101dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 102dcf2fd19SBarry Smith { 10378064530SBarry Smith PetscInt i; 10478064530SBarry Smith PetscBool identical; 10578064530SBarry Smith 106dcf2fd19SBarry Smith PetscFunctionBegin; 107dcf2fd19SBarry Smith PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1); 10878064530SBarry Smith for (i=0; i<ls->numbermonitors;i++) { 1099566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical)); 11078064530SBarry Smith if (identical) PetscFunctionReturn(0); 11178064530SBarry Smith } 1125f80ce2aSJacob Faibussowitsch PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 113dcf2fd19SBarry Smith ls->monitorftns[ls->numbermonitors] = f; 114dcf2fd19SBarry Smith ls->monitordestroy[ls->numbermonitors] = monitordestroy; 115dcf2fd19SBarry Smith ls->monitorcontext[ls->numbermonitors++] = (void*)mctx; 116dcf2fd19SBarry Smith PetscFunctionReturn(0); 117dcf2fd19SBarry Smith } 118dcf2fd19SBarry Smith 119dcf2fd19SBarry Smith /*@C 120dcf2fd19SBarry Smith SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries 121dcf2fd19SBarry Smith 122dcf2fd19SBarry Smith Collective on SNESLineSearch 123dcf2fd19SBarry Smith 124dcf2fd19SBarry Smith Input Parameters: 125dcf2fd19SBarry Smith + ls - the SNES linesearch object 126d12e167eSBarry Smith - vf - the context for the monitor, in this case it is an ASCII PetscViewer and format 127dcf2fd19SBarry Smith 128dcf2fd19SBarry Smith Level: intermediate 129dcf2fd19SBarry Smith 130db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()` 131dcf2fd19SBarry Smith @*/ 132d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf) 133dcf2fd19SBarry Smith { 134d12e167eSBarry Smith PetscViewer viewer = vf->viewer; 135dcf2fd19SBarry Smith Vec Y,W,G; 136dcf2fd19SBarry Smith 137dcf2fd19SBarry Smith PetscFunctionBegin; 1389566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,vf->format)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n")); 1419566063dSJacob Faibussowitsch PetscCall(VecView(Y,viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n")); 1439566063dSJacob Faibussowitsch PetscCall(VecView(W,viewer)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n")); 1459566063dSJacob Faibussowitsch PetscCall(VecView(G,viewer)); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 147dcf2fd19SBarry Smith PetscFunctionReturn(0); 148dcf2fd19SBarry Smith } 149dcf2fd19SBarry Smith 150f40b411bSPeter Brune /*@ 151cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 152f40b411bSPeter Brune 153cd7522eaSPeter Brune Logically Collective on Comm 154f40b411bSPeter Brune 155f40b411bSPeter Brune Input Parameters: 156cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 157f40b411bSPeter Brune 158f40b411bSPeter Brune Output Parameters: 1598e557f58SPeter Brune . outlinesearch - the new linesearch context 160f40b411bSPeter Brune 161162e0bf5SPeter Brune Level: developer 162162e0bf5SPeter Brune 163162e0bf5SPeter Brune Notes: 164162e0bf5SPeter Brune The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance 165162e0bf5SPeter Brune already associated with the SNES. This function is for developer use. 166f40b411bSPeter Brune 167db781477SPatrick Sanan .seealso: `LineSearchDestroy()`, `SNESGetLineSearch()` 168f40b411bSPeter Brune @*/ 169f40b411bSPeter Brune 170bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 171bf388a1fSBarry Smith { 172f1c6b773SPeter Brune SNESLineSearch linesearch; 173bf388a1fSBarry Smith 174bf7f4e0aSPeter Brune PetscFunctionBegin; 175ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 1769566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 1770298fd71SBarry Smith *outlinesearch = NULL; 178f5af7f23SKarl Rupp 1799566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView)); 180bf7f4e0aSPeter Brune 1810298fd71SBarry Smith linesearch->vec_sol_new = NULL; 1820298fd71SBarry Smith linesearch->vec_func_new = NULL; 1830298fd71SBarry Smith linesearch->vec_sol = NULL; 1840298fd71SBarry Smith linesearch->vec_func = NULL; 1850298fd71SBarry Smith linesearch->vec_update = NULL; 1869bd66eb0SPeter Brune 187bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 188bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 189bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 190bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 191422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 192bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 193bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 194bf7f4e0aSPeter Brune linesearch->damping = 1.0; 195bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 196bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 197516fe3c3SPeter Brune linesearch->rtol = 1e-8; 198516fe3c3SPeter Brune linesearch->atol = 1e-15; 199516fe3c3SPeter Brune linesearch->ltol = 1e-8; 2000298fd71SBarry Smith linesearch->precheckctx = NULL; 2010298fd71SBarry Smith linesearch->postcheckctx = NULL; 20259405d5eSPeter Brune linesearch->max_its = 1; 203bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 2043add74b1SBarry Smith linesearch->monitor = NULL; 205bf7f4e0aSPeter Brune *outlinesearch = linesearch; 206bf7f4e0aSPeter Brune PetscFunctionReturn(0); 207bf7f4e0aSPeter Brune } 208bf7f4e0aSPeter Brune 209f40b411bSPeter Brune /*@ 21078bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 21178bcb3b5SPeter Brune any required vectors. 212f40b411bSPeter Brune 213cd7522eaSPeter Brune Collective on SNESLineSearch 214f40b411bSPeter Brune 215f40b411bSPeter Brune Input Parameters: 216f40b411bSPeter Brune . linesearch - The LineSearch instance. 217f40b411bSPeter Brune 218cd7522eaSPeter Brune Notes: 219f190f2fcSBarry Smith For most cases, this needn't be called by users or outside of SNESLineSearchApply(). 220cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 22178bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 222cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 223cd7522eaSPeter Brune allocated upfront. 224cd7522eaSPeter Brune 22578bcb3b5SPeter Brune Level: advanced 226f40b411bSPeter Brune 227db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchReset()` 228f40b411bSPeter Brune @*/ 229f40b411bSPeter Brune 230bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 231bf388a1fSBarry Smith { 232bf7f4e0aSPeter Brune PetscFunctionBegin; 233bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 2349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 235bf7f4e0aSPeter Brune } 236bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 2379bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 2389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 2399bd66eb0SPeter Brune } 2409bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 2419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 2429bd66eb0SPeter Brune } 243*1baa6e33SBarry Smith if (linesearch->ops->setup) PetscCall((*linesearch->ops->setup)(linesearch)); 2449566063dSJacob Faibussowitsch if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch,SNESComputeFunction)); 245bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 246bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 247bf7f4e0aSPeter Brune } 248bf7f4e0aSPeter Brune PetscFunctionReturn(0); 249bf7f4e0aSPeter Brune } 250bf7f4e0aSPeter Brune 251f40b411bSPeter Brune /*@ 252f190f2fcSBarry Smith SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search. 253f40b411bSPeter Brune 254f1c6b773SPeter Brune Collective on SNESLineSearch 255f40b411bSPeter Brune 256f40b411bSPeter Brune Input Parameters: 257f40b411bSPeter Brune . linesearch - The LineSearch instance. 258f40b411bSPeter Brune 25995452b02SPatrick Sanan Notes: 26095452b02SPatrick Sanan Usually only called by SNESReset() 261f190f2fcSBarry Smith 262f190f2fcSBarry Smith Level: developer 263f40b411bSPeter Brune 264db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 265f40b411bSPeter Brune @*/ 266f40b411bSPeter Brune 267bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 268bf388a1fSBarry Smith { 269bf7f4e0aSPeter Brune PetscFunctionBegin; 2709566063dSJacob Faibussowitsch if (linesearch->ops->reset) PetscCall((*linesearch->ops->reset)(linesearch)); 271f5af7f23SKarl Rupp 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_sol_new)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&linesearch->vec_func_new)); 274bf7f4e0aSPeter Brune 2759566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 276f5af7f23SKarl Rupp 277bf7f4e0aSPeter Brune linesearch->nwork = 0; 278bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 279bf7f4e0aSPeter Brune PetscFunctionReturn(0); 280bf7f4e0aSPeter Brune } 281bf7f4e0aSPeter Brune 282ed07d7d7SPeter Brune /*@C 283f190f2fcSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search 284ed07d7d7SPeter Brune 285ed07d7d7SPeter Brune Input Parameters: 286ed07d7d7SPeter Brune . linesearch - the SNESLineSearch context 287f190f2fcSBarry Smith + func - function evaluation routine 288ed07d7d7SPeter Brune 289ed07d7d7SPeter Brune Level: developer 290ed07d7d7SPeter Brune 29195452b02SPatrick Sanan Notes: 29295452b02SPatrick Sanan This is used internally by PETSc and not called by users 293f190f2fcSBarry Smith 294db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESSetFunction()` 295ed07d7d7SPeter Brune @*/ 296ed07d7d7SPeter Brune PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec)) 297ed07d7d7SPeter Brune { 298ed07d7d7SPeter Brune PetscFunctionBegin; 299ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 300ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 301ed07d7d7SPeter Brune PetscFunctionReturn(0); 302ed07d7d7SPeter Brune } 303ed07d7d7SPeter Brune 30486d74e61SPeter Brune /*@C 305f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 306df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 307f190f2fcSBarry Smith determined the search direction. 30886d74e61SPeter Brune 309f1c6b773SPeter Brune Logically Collective on SNESLineSearch 31086d74e61SPeter Brune 31186d74e61SPeter Brune Input Parameters: 312f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 31384238204SBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence 314f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 31586d74e61SPeter Brune 31686d74e61SPeter Brune Level: intermediate 31786d74e61SPeter Brune 318db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 31986d74e61SPeter Brune @*/ 320f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx) 32186d74e61SPeter Brune { 3229bd66eb0SPeter Brune PetscFunctionBegin; 323f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 324f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 32586d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 32686d74e61SPeter Brune PetscFunctionReturn(0); 32786d74e61SPeter Brune } 32886d74e61SPeter Brune 32986d74e61SPeter Brune /*@C 33053deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 33186d74e61SPeter Brune 332f899ff85SJose E. Roman Input Parameter: 333f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 33486d74e61SPeter Brune 33586d74e61SPeter Brune Output Parameters: 33684238204SBarry Smith + func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence 337f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 33886d74e61SPeter Brune 33986d74e61SPeter Brune Level: intermediate 34086d74e61SPeter Brune 341db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 34286d74e61SPeter Brune @*/ 343f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx) 34486d74e61SPeter Brune { 3459bd66eb0SPeter Brune PetscFunctionBegin; 346f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 347f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 34886d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 34986d74e61SPeter Brune PetscFunctionReturn(0); 35086d74e61SPeter Brune } 35186d74e61SPeter Brune 35286d74e61SPeter Brune /*@C 353f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 354f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 35586d74e61SPeter Brune 356f1c6b773SPeter Brune Logically Collective on SNESLineSearch 35786d74e61SPeter Brune 35886d74e61SPeter Brune Input Parameters: 359f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 36084238204SBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPostCheck() for the calling sequence 361f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 36286d74e61SPeter Brune 36386d74e61SPeter Brune Level: intermediate 36486d74e61SPeter Brune 365db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()` 36686d74e61SPeter Brune @*/ 367f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 36886d74e61SPeter Brune { 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; 37386d74e61SPeter Brune PetscFunctionReturn(0); 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: 380f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 38186d74e61SPeter Brune 38286d74e61SPeter Brune Output Parameters: 38384238204SBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck() 384f190f2fcSBarry 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 388db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 38986d74e61SPeter Brune @*/ 390f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 39186d74e61SPeter Brune { 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; 39686d74e61SPeter Brune PetscFunctionReturn(0); 39786d74e61SPeter Brune } 39886d74e61SPeter Brune 399f40b411bSPeter Brune /*@ 400f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 401f40b411bSPeter Brune 402cd7522eaSPeter Brune Logically Collective on SNESLineSearch 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 409f40b411bSPeter Brune Output Parameters: 4108e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 411f40b411bSPeter Brune 412f190f2fcSBarry Smith Level: developer 413f40b411bSPeter Brune 414db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()` 415f40b411bSPeter Brune @*/ 4167b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 417bf7f4e0aSPeter Brune { 418bf7f4e0aSPeter Brune PetscFunctionBegin; 419bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 4206b2b7091SBarry Smith if (linesearch->ops->precheck) { 4219566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx)); 42238bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 423bf7f4e0aSPeter Brune } 424bf7f4e0aSPeter Brune PetscFunctionReturn(0); 425bf7f4e0aSPeter Brune } 426bf7f4e0aSPeter Brune 427f40b411bSPeter Brune /*@ 428ef46b1a6SStefano Zampini SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 429f40b411bSPeter Brune 430cd7522eaSPeter Brune Logically Collective on SNESLineSearch 431f40b411bSPeter Brune 432f40b411bSPeter Brune Input Parameters: 4337b1df9c1SPeter Brune + linesearch - The linesearch context 4347b1df9c1SPeter Brune . X - The last solution 4357b1df9c1SPeter Brune . Y - The step direction 4367b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 437f40b411bSPeter Brune 438f40b411bSPeter Brune Output Parameters: 43978bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 44078bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 441f40b411bSPeter Brune 442f190f2fcSBarry Smith Level: developer 443f40b411bSPeter Brune 444db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 445f40b411bSPeter Brune @*/ 4467b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 447bf7f4e0aSPeter Brune { 448bf7f4e0aSPeter Brune PetscFunctionBegin; 449bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 450bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4516b2b7091SBarry Smith if (linesearch->ops->postcheck) { 4529566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx)); 45338bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 45438bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 45586d74e61SPeter Brune } 45686d74e61SPeter Brune PetscFunctionReturn(0); 45786d74e61SPeter Brune } 45886d74e61SPeter Brune 45986d74e61SPeter Brune /*@C 46086d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 46186d74e61SPeter Brune 462cd7522eaSPeter Brune Logically Collective on SNESLineSearch 46386d74e61SPeter Brune 4644165533cSJose E. Roman Input Parameters: 46586d74e61SPeter Brune + linesearch - linesearch context 46686d74e61SPeter Brune . X - base state for this step 467907376e6SBarry Smith - ctx - context for this function 46886d74e61SPeter Brune 46997bb3fdcSJose E. Roman Input/Output Parameter: 47097bb3fdcSJose E. Roman . Y - correction, possibly modified 47197bb3fdcSJose E. Roman 47297bb3fdcSJose E. Roman Output Parameter: 47397bb3fdcSJose E. Roman . changed - flag indicating that Y was modified 47486d74e61SPeter Brune 47586d74e61SPeter Brune Options Database Key: 476cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 477cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 47886d74e61SPeter Brune 47986d74e61SPeter Brune Level: advanced 48086d74e61SPeter Brune 48186d74e61SPeter Brune Notes: 48286d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 48386d74e61SPeter Brune 48486d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 48586d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 48686d74e61SPeter Brune is generally not useful when using a Newton linearization. 48786d74e61SPeter Brune 48886d74e61SPeter Brune Reference: 48986d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 49086d74e61SPeter Brune 491db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()` 49286d74e61SPeter Brune @*/ 493f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 49486d74e61SPeter Brune { 49586d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 49686d74e61SPeter Brune Vec Ylast; 49786d74e61SPeter Brune PetscScalar dot; 49886d74e61SPeter Brune PetscInt iter; 49986d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 50086d74e61SPeter Brune SNES snes; 50186d74e61SPeter Brune 50286d74e61SPeter Brune PetscFunctionBegin; 5039566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 5049566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast)); 50586d74e61SPeter Brune if (!Ylast) { 5069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y,&Ylast)); 5079566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast)); 5089566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)Ylast)); 50986d74e61SPeter Brune } 5109566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&iter)); 51186d74e61SPeter Brune if (iter < 2) { 5129566063dSJacob Faibussowitsch PetscCall(VecCopy(Y,Ylast)); 51386d74e61SPeter Brune *changed = PETSC_FALSE; 51486d74e61SPeter Brune PetscFunctionReturn(0); 51586d74e61SPeter Brune } 51686d74e61SPeter Brune 5179566063dSJacob Faibussowitsch PetscCall(VecDot(Y,Ylast,&dot)); 5189566063dSJacob Faibussowitsch PetscCall(VecNorm(Y,NORM_2,&ynorm)); 5199566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm)); 52086d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 521255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 52286d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 52386d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 52486d74e61SPeter Brune /* Modify the step Y */ 52586d74e61SPeter Brune PetscReal alpha,ydiffnorm; 5269566063dSJacob Faibussowitsch PetscCall(VecAXPY(Ylast,-1.0,Y)); 5279566063dSJacob Faibussowitsch PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm)); 528f85e2ce2SBarry Smith alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0; 5299566063dSJacob Faibussowitsch PetscCall(VecCopy(Y,Ylast)); 5309566063dSJacob Faibussowitsch PetscCall(VecScale(Y,alpha)); 5319566063dSJacob 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)); 532a47ec85fSBarry Smith *changed = PETSC_TRUE; 53386d74e61SPeter Brune } else { 5349566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle)); 5359566063dSJacob Faibussowitsch PetscCall(VecCopy(Y,Ylast)); 53686d74e61SPeter Brune *changed = PETSC_FALSE; 537bf7f4e0aSPeter Brune } 538bf7f4e0aSPeter Brune PetscFunctionReturn(0); 539bf7f4e0aSPeter Brune } 540bf7f4e0aSPeter Brune 541f40b411bSPeter Brune /*@ 542cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 543f40b411bSPeter Brune 544f1c6b773SPeter Brune Collective on SNESLineSearch 545f40b411bSPeter Brune 546f40b411bSPeter Brune Input Parameters: 5478e557f58SPeter Brune + linesearch - The linesearch context 5488e557f58SPeter Brune - Y - The search direction 549f40b411bSPeter Brune 5506b867d5aSJose E. Roman Input/Output Parameters: 5516b867d5aSJose E. Roman + X - The current solution, on output the new solution 5526b867d5aSJose E. Roman . F - The current function, on output the new function 5536b867d5aSJose E. Roman - fnorm - The current norm, on output the new function norm 554f40b411bSPeter Brune 555cd7522eaSPeter Brune Options Database Keys: 5560b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell 557dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 5581fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 5591fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 5603c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 5613c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 562cd7522eaSPeter Brune 563cd7522eaSPeter Brune Notes: 564cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 565cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 566cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 567cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 56884238204SBarry Smith application of the line search may invoke SNESComputeFunction() several times, and 569cd7522eaSPeter Brune therefore may be fairly expensive. 570cd7522eaSPeter Brune 571f40b411bSPeter Brune Level: Intermediate 572f40b411bSPeter Brune 573db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`, 574db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetType()` 575f40b411bSPeter Brune @*/ 576bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 577bf388a1fSBarry Smith { 578bf388a1fSBarry Smith PetscFunctionBegin; 579f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 580bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 581bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 582064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 583bf7f4e0aSPeter Brune 584422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 585bf7f4e0aSPeter Brune 586bf7f4e0aSPeter Brune linesearch->vec_sol = X; 587bf7f4e0aSPeter Brune linesearch->vec_update = Y; 588bf7f4e0aSPeter Brune linesearch->vec_func = F; 589bf7f4e0aSPeter Brune 5909566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(linesearch)); 591bf7f4e0aSPeter Brune 592f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 593bf7f4e0aSPeter Brune 594f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 595f5af7f23SKarl Rupp else { 5969566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 597bf7f4e0aSPeter Brune } 598bf7f4e0aSPeter Brune 5999566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y)); 600bf7f4e0aSPeter Brune 6019566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->apply)(linesearch)); 602bf7f4e0aSPeter Brune 6039566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y)); 604bf7f4e0aSPeter Brune 605f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 606bf7f4e0aSPeter Brune PetscFunctionReturn(0); 607bf7f4e0aSPeter Brune } 608bf7f4e0aSPeter Brune 609f40b411bSPeter Brune /*@ 610f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 611f40b411bSPeter Brune 612f1c6b773SPeter Brune Collective on SNESLineSearch 613f40b411bSPeter Brune 614f40b411bSPeter Brune Input Parameters: 6158e557f58SPeter Brune . linesearch - The linesearch context 616f40b411bSPeter Brune 61784238204SBarry Smith Level: developer 618f40b411bSPeter Brune 619db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()` 620f40b411bSPeter Brune @*/ 621bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 622bf388a1fSBarry Smith { 623bf7f4e0aSPeter Brune PetscFunctionBegin; 624bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 625f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 6269e5d0892SLisandro Dalcin if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; PetscFunctionReturn(0);} 6279566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch)); 6289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset(*linesearch)); 629f5af7f23SKarl Rupp if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch); 6309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 6319566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorCancel((*linesearch))); 6329566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(linesearch)); 633bf7f4e0aSPeter Brune PetscFunctionReturn(0); 634bf7f4e0aSPeter Brune } 635bf7f4e0aSPeter Brune 636f40b411bSPeter Brune /*@ 637dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 638bf7f4e0aSPeter Brune 639bf7f4e0aSPeter Brune Input Parameters: 640dcf2fd19SBarry Smith + linesearch - the linesearch object 641dcf2fd19SBarry Smith - viewer - an ASCII PetscViewer or NULL to turn off monitor 642bf7f4e0aSPeter Brune 643dcf2fd19SBarry Smith Logically Collective on SNESLineSearch 644bf7f4e0aSPeter Brune 645bf7f4e0aSPeter Brune Options Database: 646dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor 647bf7f4e0aSPeter Brune 648bf7f4e0aSPeter Brune Level: intermediate 649bf7f4e0aSPeter Brune 650dcf2fd19SBarry Smith Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with 651d12e167eSBarry Smith SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the 652d12e167eSBarry Smith line search that are not visible to the other monitors. 653dcf2fd19SBarry Smith 654db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()` 655bf7f4e0aSPeter Brune @*/ 656dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 657bf7f4e0aSPeter Brune { 658bf7f4e0aSPeter Brune PetscFunctionBegin; 6599566063dSJacob Faibussowitsch if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer)); 6609566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&linesearch->monitor)); 661dcf2fd19SBarry Smith linesearch->monitor = viewer; 662bf7f4e0aSPeter Brune PetscFunctionReturn(0); 663bf7f4e0aSPeter Brune } 664bf7f4e0aSPeter Brune 665f40b411bSPeter Brune /*@ 666dcf2fd19SBarry Smith SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor. 6676a388c36SPeter Brune 668f190f2fcSBarry Smith Input Parameter: 6698e557f58SPeter Brune . linesearch - linesearch context 670f40b411bSPeter Brune 671f190f2fcSBarry Smith Output Parameter: 6728e557f58SPeter Brune . monitor - monitor context 673f40b411bSPeter Brune 674f40b411bSPeter Brune Logically Collective on SNES 675f40b411bSPeter Brune 676f40b411bSPeter Brune Options Database Keys: 6778e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 678f40b411bSPeter Brune 679f40b411bSPeter Brune Level: intermediate 680f40b411bSPeter Brune 681db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 682f40b411bSPeter Brune @*/ 683dcf2fd19SBarry Smith PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 6846a388c36SPeter Brune { 6856a388c36SPeter Brune PetscFunctionBegin; 686f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 6876a388c36SPeter Brune *monitor = linesearch->monitor; 6886a388c36SPeter Brune PetscFunctionReturn(0); 6896a388c36SPeter Brune } 6906a388c36SPeter Brune 691dcf2fd19SBarry Smith /*@C 692dcf2fd19SBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 693dcf2fd19SBarry Smith 694dcf2fd19SBarry Smith Collective on SNESLineSearch 695dcf2fd19SBarry Smith 696dcf2fd19SBarry Smith Input Parameters: 697dcf2fd19SBarry Smith + ls - LineSearch object you wish to monitor 698dcf2fd19SBarry Smith . name - the monitor type one is seeking 699dcf2fd19SBarry Smith . help - message indicating what monitoring is done 700dcf2fd19SBarry Smith . manual - manual page for the monitor 701dcf2fd19SBarry Smith . monitor - the monitor function 702dcf2fd19SBarry 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 objects 703dcf2fd19SBarry Smith 704dcf2fd19SBarry Smith Level: developer 705dcf2fd19SBarry Smith 706db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 707db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 708db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 709db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 710c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 711db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 712db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 713dcf2fd19SBarry Smith @*/ 714d12e167eSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*)) 715dcf2fd19SBarry Smith { 716dcf2fd19SBarry Smith PetscViewer viewer; 717dcf2fd19SBarry Smith PetscViewerFormat format; 718dcf2fd19SBarry Smith PetscBool flg; 719dcf2fd19SBarry Smith 720dcf2fd19SBarry Smith PetscFunctionBegin; 7219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg)); 722dcf2fd19SBarry Smith if (flg) { 723d12e167eSBarry Smith PetscViewerAndFormat *vf; 7249566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf)); 7259566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 726*1baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls,vf)); 7279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy)); 728dcf2fd19SBarry Smith } 729dcf2fd19SBarry Smith PetscFunctionReturn(0); 730dcf2fd19SBarry Smith } 731dcf2fd19SBarry Smith 732f40b411bSPeter Brune /*@ 733f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 734f40b411bSPeter Brune 735f899ff85SJose E. Roman Input Parameter: 7368e557f58SPeter Brune . linesearch - linesearch context 737f40b411bSPeter Brune 738f40b411bSPeter Brune Options Database Keys: 7390b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 7403c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 7411fe24845SBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms()) 74271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 7431a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 7441a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 7451a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 7461a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 7471a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 748dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches 749dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 7508e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 751cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 7521a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 753d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 754f40b411bSPeter Brune 755f1c6b773SPeter Brune Logically Collective on SNESLineSearch 756f40b411bSPeter Brune 757f40b411bSPeter Brune Level: intermediate 758f40b411bSPeter Brune 759db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 760db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 761f40b411bSPeter Brune @*/ 762bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 763bf388a1fSBarry Smith { 7641a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 765bf7f4e0aSPeter Brune char type[256]; 766bf7f4e0aSPeter Brune PetscBool flg, set; 767dcf2fd19SBarry Smith PetscViewer viewer; 768bf388a1fSBarry Smith 769bf7f4e0aSPeter Brune PetscFunctionBegin; 7709566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll()); 771bf7f4e0aSPeter Brune 772d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch); 773f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 7749566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg)); 775bf7f4e0aSPeter Brune if (flg) { 7769566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,type)); 777bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 7789566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,deft)); 779bf7f4e0aSPeter Brune } 780bf7f4e0aSPeter Brune 7819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set)); 782dcf2fd19SBarry Smith if (set) { 7839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetDefaultMonitor(linesearch,viewer)); 7849566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 785dcf2fd19SBarry Smith } 7869566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL)); 787bf7f4e0aSPeter Brune 7881a9b3a06SPeter Brune /* tolerances */ 7899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL)); 7909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL)); 7919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL)); 7929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL)); 7939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL)); 7949566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL)); 7957a35526eSPeter Brune 7961a9b3a06SPeter Brune /* damping parameters */ 7979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL)); 7981a9b3a06SPeter Brune 7999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL)); 8001a9b3a06SPeter Brune 8011a9b3a06SPeter Brune /* precheck */ 8029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set)); 8037a35526eSPeter Brune if (set) { 8047a35526eSPeter Brune if (flg) { 8057a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 806f5af7f23SKarl Rupp 807d0609cedSBarry 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)); 8089566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle)); 8097a35526eSPeter Brune } else { 8109566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch,NULL,NULL)); 8117a35526eSPeter Brune } 8127a35526eSPeter Brune } 8139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL)); 8149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL)); 8157a35526eSPeter Brune 816*1baa6e33SBarry Smith if (linesearch->ops->setfromoptions) PetscCall((*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch)); 8175a9b6599SPeter Brune 8189566063dSJacob Faibussowitsch PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch)); 819d0609cedSBarry Smith PetscOptionsEnd(); 820bf7f4e0aSPeter Brune PetscFunctionReturn(0); 821bf7f4e0aSPeter Brune } 822bf7f4e0aSPeter Brune 823f40b411bSPeter Brune /*@ 824f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 825f40b411bSPeter Brune 826f40b411bSPeter Brune Input Parameters: 8278e557f58SPeter Brune . linesearch - linesearch context 828f40b411bSPeter Brune 829f1c6b773SPeter Brune Logically Collective on SNESLineSearch 830f40b411bSPeter Brune 831f40b411bSPeter Brune Level: intermediate 832f40b411bSPeter Brune 833db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()` 834f40b411bSPeter Brune @*/ 835bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 836bf388a1fSBarry Smith { 8377f1410a3SPeter Brune PetscBool iascii; 838bf388a1fSBarry Smith 839bf7f4e0aSPeter Brune PetscFunctionBegin; 8407f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 8417f1410a3SPeter Brune if (!viewer) { 8429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer)); 8437f1410a3SPeter Brune } 8447f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 8457f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 846f40b411bSPeter Brune 8479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 8487f1410a3SPeter Brune if (iascii) { 8499566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer)); 8507f1410a3SPeter Brune if (linesearch->ops->view) { 8519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8529566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->view)(linesearch,viewer)); 8539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 8547f1410a3SPeter Brune } 8559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol)); 8569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol)); 85763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 8586b2b7091SBarry Smith if (linesearch->ops->precheck) { 8596b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 86063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n")); 8617f1410a3SPeter Brune } else { 86263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n")); 8637f1410a3SPeter Brune } 8647f1410a3SPeter Brune } 8656b2b7091SBarry Smith if (linesearch->ops->postcheck) { 86663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n")); 8677f1410a3SPeter Brune } 8687f1410a3SPeter Brune } 869bf7f4e0aSPeter Brune PetscFunctionReturn(0); 870bf7f4e0aSPeter Brune } 871bf7f4e0aSPeter Brune 872ea5d4fccSPeter Brune /*@C 873a80ff896SJed Brown SNESLineSearchGetType - Gets the linesearch type 874a80ff896SJed Brown 875a80ff896SJed Brown Logically Collective on SNESLineSearch 876a80ff896SJed Brown 877a80ff896SJed Brown Input Parameters: 878a80ff896SJed Brown . linesearch - linesearch context 879a80ff896SJed Brown 880a80ff896SJed Brown Output Parameters: 881a80ff896SJed Brown - type - The type of line search, or NULL if not set 882a80ff896SJed Brown 883a80ff896SJed Brown Level: intermediate 884a80ff896SJed Brown 885db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 886a80ff896SJed Brown @*/ 887a80ff896SJed Brown PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 888a80ff896SJed Brown { 889a80ff896SJed Brown PetscFunctionBegin; 890a80ff896SJed Brown PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 891dadcf809SJacob Faibussowitsch PetscValidPointer(type,2); 892a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name; 893a80ff896SJed Brown PetscFunctionReturn(0); 894a80ff896SJed Brown } 895a80ff896SJed Brown 896a80ff896SJed Brown /*@C 897f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 898f40b411bSPeter Brune 899f190f2fcSBarry Smith Logically Collective on SNESLineSearch 900f190f2fcSBarry Smith 901f40b411bSPeter Brune Input Parameters: 9028e557f58SPeter Brune + linesearch - linesearch context 903f40b411bSPeter Brune - type - The type of line search to be used 904f40b411bSPeter Brune 905cd7522eaSPeter Brune Available Types: 9060b00b554SBarry Smith + SNESLINESEARCHBASIC - (or equivalently SNESLINESEARCHNONE) Simple damping line search, defaults to using the full Newton step 9071fe24845SBarry Smith . SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function 9081fe24845SBarry Smith . SNESLINESEARCHL2 - Secant line search over the L2 norm of the function 9091fe24845SBarry Smith . SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 9101fe24845SBarry Smith . SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch 9111fe24845SBarry Smith - SNESLINESEARCHSHELL - User provided SNESLineSearch implementation 9121fe24845SBarry Smith 9131fe24845SBarry Smith Options Database: 9140b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 915cd7522eaSPeter Brune 916f40b411bSPeter Brune Level: intermediate 917f40b411bSPeter Brune 918db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 919f40b411bSPeter Brune @*/ 92019fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 921bf7f4e0aSPeter Brune { 922bf7f4e0aSPeter Brune PetscBool match; 9235f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch); 924bf7f4e0aSPeter Brune 925bf7f4e0aSPeter Brune PetscFunctionBegin; 926f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 927bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 928bf7f4e0aSPeter Brune 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch,type,&match)); 930bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 931bf7f4e0aSPeter Brune 9329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList,type,&r)); 9335f80ce2aSJacob Faibussowitsch PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 934bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 935bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 9369566063dSJacob Faibussowitsch PetscCall((*(linesearch)->ops->destroy)(linesearch)); 9370298fd71SBarry Smith linesearch->ops->destroy = NULL; 938bf7f4e0aSPeter Brune } 939f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 9409e5d0892SLisandro Dalcin linesearch->ops->apply = NULL; 9419e5d0892SLisandro Dalcin linesearch->ops->view = NULL; 9429e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL; 9439e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL; 944bf7f4e0aSPeter Brune 9459566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch,type)); 9469566063dSJacob Faibussowitsch PetscCall((*r)(linesearch)); 947bf7f4e0aSPeter Brune PetscFunctionReturn(0); 948bf7f4e0aSPeter Brune } 949bf7f4e0aSPeter Brune 950f40b411bSPeter Brune /*@ 95178bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 952f40b411bSPeter Brune 953f40b411bSPeter Brune Input Parameters: 9548e557f58SPeter Brune + linesearch - linesearch context 955f40b411bSPeter Brune - snes - The snes instance 956f40b411bSPeter Brune 95778bcb3b5SPeter Brune Level: developer 95878bcb3b5SPeter Brune 95978bcb3b5SPeter Brune Notes: 960f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 9617601faf0SJed Brown SNESGetLineSearch(). This routine is therefore mainly called within SNES 96278bcb3b5SPeter Brune implementations. 963f40b411bSPeter Brune 9648141a3b9SPeter Brune Level: developer 965f40b411bSPeter Brune 966db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 967f40b411bSPeter Brune @*/ 968bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 969bf388a1fSBarry Smith { 970bf7f4e0aSPeter Brune PetscFunctionBegin; 971f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 972bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 973bf7f4e0aSPeter Brune linesearch->snes = snes; 974bf7f4e0aSPeter Brune PetscFunctionReturn(0); 975bf7f4e0aSPeter Brune } 976bf7f4e0aSPeter Brune 977f40b411bSPeter Brune /*@ 9788141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 9798141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 9808141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 9818141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 982f40b411bSPeter Brune 983f40b411bSPeter Brune Input Parameters: 9848e557f58SPeter Brune . linesearch - linesearch context 985f40b411bSPeter Brune 986f40b411bSPeter Brune Output Parameters: 987f40b411bSPeter Brune . snes - The snes instance 988f40b411bSPeter Brune 9898141a3b9SPeter Brune Level: developer 990f40b411bSPeter Brune 991db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 992f40b411bSPeter Brune @*/ 993bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 994bf388a1fSBarry Smith { 995bf7f4e0aSPeter Brune PetscFunctionBegin; 996f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9976a388c36SPeter Brune PetscValidPointer(snes,2); 998bf7f4e0aSPeter Brune *snes = linesearch->snes; 999bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1000bf7f4e0aSPeter Brune } 1001bf7f4e0aSPeter Brune 1002f40b411bSPeter Brune /*@ 1003f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 1004f40b411bSPeter Brune 1005f40b411bSPeter Brune Input Parameters: 10068e557f58SPeter Brune . linesearch - linesearch context 1007f40b411bSPeter Brune 1008f40b411bSPeter Brune Output Parameters: 1009cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 1010f40b411bSPeter Brune 101178bcb3b5SPeter Brune Level: advanced 101278bcb3b5SPeter Brune 10138e557f58SPeter Brune Notes: 10148e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 101578bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 101678bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 101778bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 101878bcb3b5SPeter Brune 1019db781477SPatrick Sanan .seealso: `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1020f40b411bSPeter Brune @*/ 1021f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 10226a388c36SPeter Brune { 10236a388c36SPeter Brune PetscFunctionBegin; 1024f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1025534a8f05SLisandro Dalcin PetscValidRealPointer(lambda, 2); 10266a388c36SPeter Brune *lambda = linesearch->lambda; 10276a388c36SPeter Brune PetscFunctionReturn(0); 10286a388c36SPeter Brune } 10296a388c36SPeter Brune 1030f40b411bSPeter Brune /*@ 1031f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 1032f40b411bSPeter Brune 1033f40b411bSPeter Brune Input Parameters: 10348e557f58SPeter Brune + linesearch - linesearch context 1035f40b411bSPeter Brune - lambda - The last steplength. 1036f40b411bSPeter Brune 1037cd7522eaSPeter Brune Notes: 1038f190f2fcSBarry Smith This routine is typically used within implementations of SNESLineSearchApply() 1039cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 1040cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 1041cd7522eaSPeter Brune as an inner scaling parameter. 1042cd7522eaSPeter Brune 104378bcb3b5SPeter Brune Level: advanced 1044f40b411bSPeter Brune 1045db781477SPatrick Sanan .seealso: `SNESLineSearchGetLambda()` 1046f40b411bSPeter Brune @*/ 1047f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 10486a388c36SPeter Brune { 10496a388c36SPeter Brune PetscFunctionBegin; 1050f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 10516a388c36SPeter Brune linesearch->lambda = lambda; 10526a388c36SPeter Brune PetscFunctionReturn(0); 10536a388c36SPeter Brune } 10546a388c36SPeter Brune 1055f40b411bSPeter Brune /*@ 10563c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 105778bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 105878bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 105978bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1060f40b411bSPeter Brune 1061f899ff85SJose E. Roman Input Parameter: 10628e557f58SPeter Brune . linesearch - linesearch context 1063f40b411bSPeter Brune 1064f40b411bSPeter Brune Output Parameters: 1065516fe3c3SPeter Brune + steptol - The minimum steplength 10666cc8e53bSPeter Brune . maxstep - The maximum steplength 1067516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1068516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1069516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1070516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1071f40b411bSPeter Brune 107278bcb3b5SPeter Brune Level: intermediate 107378bcb3b5SPeter Brune 107478bcb3b5SPeter Brune Notes: 107578bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 10763c7d6663SPeter Brune the type requires. 1077516fe3c3SPeter Brune 1078db781477SPatrick Sanan .seealso: `SNESLineSearchSetTolerances()` 1079f40b411bSPeter Brune @*/ 1080f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 10816a388c36SPeter Brune { 10826a388c36SPeter Brune PetscFunctionBegin; 1083f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1084516fe3c3SPeter Brune if (steptol) { 1085534a8f05SLisandro Dalcin PetscValidRealPointer(steptol, 2); 10866a388c36SPeter Brune *steptol = linesearch->steptol; 1087516fe3c3SPeter Brune } 1088516fe3c3SPeter Brune if (maxstep) { 1089534a8f05SLisandro Dalcin PetscValidRealPointer(maxstep, 3); 1090516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1091516fe3c3SPeter Brune } 1092516fe3c3SPeter Brune if (rtol) { 1093534a8f05SLisandro Dalcin PetscValidRealPointer(rtol, 4); 1094516fe3c3SPeter Brune *rtol = linesearch->rtol; 1095516fe3c3SPeter Brune } 1096516fe3c3SPeter Brune if (atol) { 1097534a8f05SLisandro Dalcin PetscValidRealPointer(atol, 5); 1098516fe3c3SPeter Brune *atol = linesearch->atol; 1099516fe3c3SPeter Brune } 1100516fe3c3SPeter Brune if (ltol) { 1101534a8f05SLisandro Dalcin PetscValidRealPointer(ltol, 6); 1102516fe3c3SPeter Brune *ltol = linesearch->ltol; 1103516fe3c3SPeter Brune } 1104516fe3c3SPeter Brune if (max_its) { 1105534a8f05SLisandro Dalcin PetscValidIntPointer(max_its, 7); 1106516fe3c3SPeter Brune *max_its = linesearch->max_its; 1107516fe3c3SPeter Brune } 11086a388c36SPeter Brune PetscFunctionReturn(0); 11096a388c36SPeter Brune } 11106a388c36SPeter Brune 1111f40b411bSPeter Brune /*@ 11123c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 111378bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 111478bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 111578bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1116f40b411bSPeter Brune 1117f40b411bSPeter Brune Input Parameters: 11188e557f58SPeter Brune + linesearch - linesearch context 1119516fe3c3SPeter Brune . steptol - The minimum steplength 11206cc8e53bSPeter Brune . maxstep - The maximum steplength 1121516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1122516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1123516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1124516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1125f40b411bSPeter Brune 112678bcb3b5SPeter Brune Notes: 11273c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 112878bcb3b5SPeter Brune place of an argument. 1129f40b411bSPeter Brune 113078bcb3b5SPeter Brune Level: intermediate 1131516fe3c3SPeter Brune 1132db781477SPatrick Sanan .seealso: `SNESLineSearchGetTolerances()` 1133f40b411bSPeter Brune @*/ 1134f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 11356a388c36SPeter Brune { 11366a388c36SPeter Brune PetscFunctionBegin; 1137f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1138d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1139d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1140d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1141d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1142d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1143d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1144d3952184SSatish Balay 1145d3952184SSatish Balay if (steptol!= PETSC_DEFAULT) { 11465f80ce2aSJacob Faibussowitsch PetscCheck(steptol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 11476a388c36SPeter Brune linesearch->steptol = steptol; 1148d3952184SSatish Balay } 1149d3952184SSatish Balay 1150d3952184SSatish Balay if (maxstep!= PETSC_DEFAULT) { 11515f80ce2aSJacob Faibussowitsch PetscCheck(maxstep >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1152516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1153d3952184SSatish Balay } 1154d3952184SSatish Balay 1155d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 11562061ca32SJunchao 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); 1157516fe3c3SPeter Brune linesearch->rtol = rtol; 1158d3952184SSatish Balay } 1159d3952184SSatish Balay 1160d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 11615f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1162516fe3c3SPeter Brune linesearch->atol = atol; 1163d3952184SSatish Balay } 1164d3952184SSatish Balay 1165d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 11665f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol); 1167516fe3c3SPeter Brune linesearch->ltol = ltol; 1168d3952184SSatish Balay } 1169d3952184SSatish Balay 1170d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 117163a3b9bcSJacob Faibussowitsch PetscCheck(max_its >= 0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %" PetscInt_FMT " must be non-negative",max_its); 1172516fe3c3SPeter Brune linesearch->max_its = max_its; 1173d3952184SSatish Balay } 11746a388c36SPeter Brune PetscFunctionReturn(0); 11756a388c36SPeter Brune } 11766a388c36SPeter Brune 1177f40b411bSPeter Brune /*@ 1178f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1179f40b411bSPeter Brune 1180f40b411bSPeter Brune Input Parameters: 11818e557f58SPeter Brune . linesearch - linesearch context 1182f40b411bSPeter Brune 1183f40b411bSPeter Brune Output Parameters: 11848e557f58SPeter Brune . damping - The damping parameter 1185f40b411bSPeter Brune 118678bcb3b5SPeter Brune Level: advanced 1187f40b411bSPeter Brune 1188db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN` 1189f40b411bSPeter Brune @*/ 1190f40b411bSPeter Brune 1191f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 11926a388c36SPeter Brune { 11936a388c36SPeter Brune PetscFunctionBegin; 1194f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1195534a8f05SLisandro Dalcin PetscValidRealPointer(damping, 2); 11966a388c36SPeter Brune *damping = linesearch->damping; 11976a388c36SPeter Brune PetscFunctionReturn(0); 11986a388c36SPeter Brune } 11996a388c36SPeter Brune 1200f40b411bSPeter Brune /*@ 1201fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter. 1202f40b411bSPeter Brune 1203f40b411bSPeter Brune Input Parameters: 120403fd4120SBarry Smith + linesearch - linesearch context 120503fd4120SBarry Smith - damping - The damping parameter 1206f40b411bSPeter Brune 120703fd4120SBarry Smith Options Database: 120803fd4120SBarry Smith . -snes_linesearch_damping 1209f40b411bSPeter Brune Level: intermediate 1210f40b411bSPeter Brune 1211cd7522eaSPeter Brune Notes: 12120b00b554SBarry Smith The basic (also known as the none) line search merely takes the update step scaled by the damping parameter. 1213cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 121478bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1215cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1216cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1217cd7522eaSPeter Brune 1218db781477SPatrick Sanan .seealso: `SNESLineSearchGetDamping()` 1219f40b411bSPeter Brune @*/ 1220f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 12216a388c36SPeter Brune { 12226a388c36SPeter Brune PetscFunctionBegin; 1223f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12246a388c36SPeter Brune linesearch->damping = damping; 12256a388c36SPeter Brune PetscFunctionReturn(0); 12266a388c36SPeter Brune } 12276a388c36SPeter Brune 122859405d5eSPeter Brune /*@ 122959405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 123059405d5eSPeter Brune 123159405d5eSPeter Brune Input Parameters: 123278bcb3b5SPeter Brune . linesearch - linesearch context 123359405d5eSPeter Brune 123459405d5eSPeter Brune Output Parameters: 12358e557f58SPeter Brune . order - The order 123659405d5eSPeter Brune 123778bcb3b5SPeter Brune Possible Values for order: 12383c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 12393c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 12403c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 124178bcb3b5SPeter Brune 124259405d5eSPeter Brune Level: intermediate 124359405d5eSPeter Brune 1244db781477SPatrick Sanan .seealso: `SNESLineSearchSetOrder()` 124559405d5eSPeter Brune @*/ 124659405d5eSPeter Brune 1247b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 124859405d5eSPeter Brune { 124959405d5eSPeter Brune PetscFunctionBegin; 125059405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1251534a8f05SLisandro Dalcin PetscValidIntPointer(order, 2); 125259405d5eSPeter Brune *order = linesearch->order; 125359405d5eSPeter Brune PetscFunctionReturn(0); 125459405d5eSPeter Brune } 125559405d5eSPeter Brune 125659405d5eSPeter Brune /*@ 12571f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 125859405d5eSPeter Brune 125959405d5eSPeter Brune Input Parameters: 1260a2b725a8SWilliam Gropp + linesearch - linesearch context 1261a2b725a8SWilliam Gropp - order - The damping parameter 126259405d5eSPeter Brune 126359405d5eSPeter Brune Level: intermediate 126459405d5eSPeter Brune 126578bcb3b5SPeter Brune Possible Values for order: 12663c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 12673c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 12683c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 126978bcb3b5SPeter Brune 127059405d5eSPeter Brune Notes: 127159405d5eSPeter Brune Variable orders are supported by the following line searches: 127278bcb3b5SPeter Brune + bt - cubic and quadratic 127378bcb3b5SPeter Brune - cp - linear and quadratic 127459405d5eSPeter Brune 1275db781477SPatrick Sanan .seealso: `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 127659405d5eSPeter Brune @*/ 1277b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 127859405d5eSPeter Brune { 127959405d5eSPeter Brune PetscFunctionBegin; 128059405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 128159405d5eSPeter Brune linesearch->order = order; 128259405d5eSPeter Brune PetscFunctionReturn(0); 128359405d5eSPeter Brune } 128459405d5eSPeter Brune 1285f40b411bSPeter Brune /*@ 1286f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1287f40b411bSPeter Brune 1288f899ff85SJose E. Roman Input Parameter: 128978bcb3b5SPeter Brune . linesearch - linesearch context 1290f40b411bSPeter Brune 1291f40b411bSPeter Brune Output Parameters: 1292f40b411bSPeter Brune + xnorm - The norm of the current solution 1293f40b411bSPeter Brune . fnorm - The norm of the current function 1294f40b411bSPeter Brune - ynorm - The norm of the current update 1295f40b411bSPeter Brune 1296cd7522eaSPeter Brune Notes: 1297cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1298cd7522eaSPeter Brune 129978bcb3b5SPeter Brune Level: developer 1300f40b411bSPeter Brune 1301db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1302f40b411bSPeter Brune @*/ 1303f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1304bf7f4e0aSPeter Brune { 1305bf7f4e0aSPeter Brune PetscFunctionBegin; 1306f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1307f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1308f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1309f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 1310bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1311bf7f4e0aSPeter Brune } 1312bf7f4e0aSPeter Brune 1313f40b411bSPeter Brune /*@ 1314f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1315f40b411bSPeter Brune 1316f40b411bSPeter Brune Input Parameters: 131778bcb3b5SPeter Brune + linesearch - linesearch context 1318f40b411bSPeter Brune . xnorm - The norm of the current solution 1319f40b411bSPeter Brune . fnorm - The norm of the current function 1320f40b411bSPeter Brune - ynorm - The norm of the current update 1321f40b411bSPeter Brune 132278bcb3b5SPeter Brune Level: advanced 1323f40b411bSPeter Brune 1324db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1325f40b411bSPeter Brune @*/ 1326f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 13276a388c36SPeter Brune { 13286a388c36SPeter Brune PetscFunctionBegin; 1329f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13306a388c36SPeter Brune linesearch->xnorm = xnorm; 13316a388c36SPeter Brune linesearch->fnorm = fnorm; 13326a388c36SPeter Brune linesearch->ynorm = ynorm; 13336a388c36SPeter Brune PetscFunctionReturn(0); 13346a388c36SPeter Brune } 13356a388c36SPeter Brune 1336f40b411bSPeter Brune /*@ 1337f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1338f40b411bSPeter Brune 1339f40b411bSPeter Brune Input Parameters: 134078bcb3b5SPeter Brune . linesearch - linesearch context 1341f40b411bSPeter Brune 1342f40b411bSPeter Brune Options Database Keys: 13438e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1344f40b411bSPeter Brune 1345f40b411bSPeter Brune Level: intermediate 1346f40b411bSPeter Brune 1347db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1348f40b411bSPeter Brune @*/ 1349f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 13506a388c36SPeter Brune { 13519bd66eb0SPeter Brune SNES snes; 1352bf388a1fSBarry Smith 13536a388c36SPeter Brune PetscFunctionBegin; 13546a388c36SPeter Brune if (linesearch->norms) { 13559bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 13569566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 13579566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13589566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13599566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 13609bd66eb0SPeter Brune } else { 13619566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13629566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13639566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13649566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 13659566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 13669566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 13676a388c36SPeter Brune } 13689bd66eb0SPeter Brune } 13696a388c36SPeter Brune PetscFunctionReturn(0); 13706a388c36SPeter Brune } 13716a388c36SPeter Brune 13726f263ca3SPeter Brune /*@ 13736f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 13746f263ca3SPeter Brune 13756f263ca3SPeter Brune Input Parameters: 137678bcb3b5SPeter Brune + linesearch - linesearch context 137778bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 13786f263ca3SPeter Brune 13796f263ca3SPeter Brune Options Database Keys: 13800b00b554SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch 13816f263ca3SPeter Brune 13826f263ca3SPeter Brune Notes: 13830b00b554SBarry 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. 13846f263ca3SPeter Brune 13856f263ca3SPeter Brune Level: intermediate 13866f263ca3SPeter Brune 1387db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 13886f263ca3SPeter Brune @*/ 13896f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 13906f263ca3SPeter Brune { 13916f263ca3SPeter Brune PetscFunctionBegin; 13926f263ca3SPeter Brune linesearch->norms = flg; 13936f263ca3SPeter Brune PetscFunctionReturn(0); 13946f263ca3SPeter Brune } 13956f263ca3SPeter Brune 1396f40b411bSPeter Brune /*@ 1397f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1398f40b411bSPeter Brune 1399f899ff85SJose E. Roman Input Parameter: 140078bcb3b5SPeter Brune . linesearch - linesearch context 1401f40b411bSPeter Brune 1402f40b411bSPeter Brune Output Parameters: 14036232e825SPeter Brune + X - Solution vector 14046232e825SPeter Brune . F - Function vector 14056232e825SPeter Brune . Y - Search direction vector 14066232e825SPeter Brune . W - Solution work vector 14076232e825SPeter Brune - G - Function work vector 14086232e825SPeter Brune 14097bba9028SPeter Brune Notes: 14107bba9028SPeter Brune At the beginning of a line search application, X should contain a 14116232e825SPeter Brune solution and the vector F the function computed at X. At the end of the 14126232e825SPeter Brune line search application, X should contain the new solution, and F the 14136232e825SPeter Brune function evaluated at the new solution. 1414f40b411bSPeter Brune 14152a7a6963SBarry Smith These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 14162a7a6963SBarry Smith 141778bcb3b5SPeter Brune Level: advanced 1418f40b411bSPeter Brune 1419db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1420f40b411bSPeter Brune @*/ 1421bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1422bf388a1fSBarry Smith { 14236a388c36SPeter Brune PetscFunctionBegin; 1424f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14256a388c36SPeter Brune if (X) { 14266a388c36SPeter Brune PetscValidPointer(X, 2); 14276a388c36SPeter Brune *X = linesearch->vec_sol; 14286a388c36SPeter Brune } 14296a388c36SPeter Brune if (F) { 14306a388c36SPeter Brune PetscValidPointer(F, 3); 14316a388c36SPeter Brune *F = linesearch->vec_func; 14326a388c36SPeter Brune } 14336a388c36SPeter Brune if (Y) { 14346a388c36SPeter Brune PetscValidPointer(Y, 4); 14356a388c36SPeter Brune *Y = linesearch->vec_update; 14366a388c36SPeter Brune } 14376a388c36SPeter Brune if (W) { 14386a388c36SPeter Brune PetscValidPointer(W, 5); 14396a388c36SPeter Brune *W = linesearch->vec_sol_new; 14406a388c36SPeter Brune } 14416a388c36SPeter Brune if (G) { 14426a388c36SPeter Brune PetscValidPointer(G, 6); 14436a388c36SPeter Brune *G = linesearch->vec_func_new; 14446a388c36SPeter Brune } 14456a388c36SPeter Brune PetscFunctionReturn(0); 14466a388c36SPeter Brune } 14476a388c36SPeter Brune 1448f40b411bSPeter Brune /*@ 1449f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1450f40b411bSPeter Brune 1451f40b411bSPeter Brune Input Parameters: 145278bcb3b5SPeter Brune + linesearch - linesearch context 14536232e825SPeter Brune . X - Solution vector 14546232e825SPeter Brune . F - Function vector 14556232e825SPeter Brune . Y - Search direction vector 14566232e825SPeter Brune . W - Solution work vector 14576232e825SPeter Brune - G - Function work vector 1458f40b411bSPeter Brune 145978bcb3b5SPeter Brune Level: advanced 1460f40b411bSPeter Brune 1461db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1462f40b411bSPeter Brune @*/ 1463bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1464bf388a1fSBarry Smith { 14656a388c36SPeter Brune PetscFunctionBegin; 1466f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14676a388c36SPeter Brune if (X) { 14686a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 14696a388c36SPeter Brune linesearch->vec_sol = X; 14706a388c36SPeter Brune } 14716a388c36SPeter Brune if (F) { 14726a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 14736a388c36SPeter Brune linesearch->vec_func = F; 14746a388c36SPeter Brune } 14756a388c36SPeter Brune if (Y) { 14766a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 14776a388c36SPeter Brune linesearch->vec_update = Y; 14786a388c36SPeter Brune } 14796a388c36SPeter Brune if (W) { 14806a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 14816a388c36SPeter Brune linesearch->vec_sol_new = W; 14826a388c36SPeter Brune } 14836a388c36SPeter Brune if (G) { 14846a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 14856a388c36SPeter Brune linesearch->vec_func_new = G; 14866a388c36SPeter Brune } 14876a388c36SPeter Brune PetscFunctionReturn(0); 14886a388c36SPeter Brune } 14896a388c36SPeter Brune 1490e7058c64SPeter Brune /*@C 1491f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1492e7058c64SPeter Brune SNES options in the database. 1493e7058c64SPeter Brune 1494cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1495e7058c64SPeter Brune 1496e7058c64SPeter Brune Input Parameters: 1497e7058c64SPeter Brune + snes - the SNES context 1498e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1499e7058c64SPeter Brune 1500e7058c64SPeter Brune Notes: 1501e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1502e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1503e7058c64SPeter Brune 1504e7058c64SPeter Brune Level: advanced 1505e7058c64SPeter Brune 1506db781477SPatrick Sanan .seealso: `SNESGetOptionsPrefix()` 1507e7058c64SPeter Brune @*/ 1508f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1509e7058c64SPeter Brune { 1510e7058c64SPeter Brune PetscFunctionBegin; 1511f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15129566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix)); 1513e7058c64SPeter Brune PetscFunctionReturn(0); 1514e7058c64SPeter Brune } 1515e7058c64SPeter Brune 1516e7058c64SPeter Brune /*@C 1517f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1518f1c6b773SPeter Brune SNESLineSearch options in the database. 1519e7058c64SPeter Brune 1520e7058c64SPeter Brune Not Collective 1521e7058c64SPeter Brune 1522e7058c64SPeter Brune Input Parameter: 1523cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1524e7058c64SPeter Brune 1525e7058c64SPeter Brune Output Parameter: 1526e7058c64SPeter Brune . prefix - pointer to the prefix string used 1527e7058c64SPeter Brune 15288e557f58SPeter Brune Notes: 15298e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1530e7058c64SPeter Brune sufficient length to hold the prefix. 1531e7058c64SPeter Brune 1532e7058c64SPeter Brune Level: advanced 1533e7058c64SPeter Brune 1534db781477SPatrick Sanan .seealso: `SNESAppendOptionsPrefix()` 1535e7058c64SPeter Brune @*/ 1536f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1537e7058c64SPeter Brune { 1538e7058c64SPeter Brune PetscFunctionBegin; 1539f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 15409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix)); 1541e7058c64SPeter Brune PetscFunctionReturn(0); 1542e7058c64SPeter Brune } 1543bf7f4e0aSPeter Brune 15448d359177SBarry Smith /*@C 1545fa0ddf94SBarry Smith SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1546f40b411bSPeter Brune 1547d8d19677SJose E. Roman Input Parameters: 1548f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1549f40b411bSPeter Brune - nwork - the number of work vectors 1550f40b411bSPeter Brune 1551f40b411bSPeter Brune Level: developer 1552f40b411bSPeter Brune 1553db781477SPatrick Sanan .seealso: `SNESSetWorkVecs()` 1554f40b411bSPeter Brune @*/ 1555fa0ddf94SBarry Smith PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1556bf7f4e0aSPeter Brune { 1557bf7f4e0aSPeter Brune PetscFunctionBegin; 1558bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 15599566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 15608d359177SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1561bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1562bf7f4e0aSPeter Brune } 1563bf7f4e0aSPeter Brune 1564f40b411bSPeter Brune /*@ 1565422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1566f40b411bSPeter Brune 1567f40b411bSPeter Brune Input Parameters: 156878bcb3b5SPeter Brune . linesearch - linesearch context 1569f40b411bSPeter Brune 1570f40b411bSPeter Brune Output Parameters: 1571422a814eSBarry Smith . result - The success or failure status 1572f40b411bSPeter Brune 1573cd7522eaSPeter Brune Notes: 1574c98378a5SBarry Smith This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1575cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1576cd7522eaSPeter Brune 1577f40b411bSPeter Brune Level: intermediate 1578f40b411bSPeter Brune 1579db781477SPatrick Sanan .seealso: `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1580f40b411bSPeter Brune @*/ 1581422a814eSBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1582bf7f4e0aSPeter Brune { 1583bf7f4e0aSPeter Brune PetscFunctionBegin; 1584f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1585422a814eSBarry Smith PetscValidPointer(result, 2); 1586422a814eSBarry Smith *result = linesearch->result; 1587bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1588bf7f4e0aSPeter Brune } 1589bf7f4e0aSPeter Brune 1590f40b411bSPeter Brune /*@ 1591422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1592f40b411bSPeter Brune 1593f40b411bSPeter Brune Input Parameters: 159478bcb3b5SPeter Brune + linesearch - linesearch context 1595422a814eSBarry Smith - result - The success or failure status 1596f40b411bSPeter Brune 1597cd7522eaSPeter Brune Notes: 1598cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1599cd7522eaSPeter Brune the success or failure of the line search method. 1600cd7522eaSPeter Brune 160178bcb3b5SPeter Brune Level: developer 1602f40b411bSPeter Brune 1603db781477SPatrick Sanan .seealso: `SNESLineSearchGetSResult()` 1604f40b411bSPeter Brune @*/ 1605422a814eSBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 16066a388c36SPeter Brune { 16076a388c36SPeter Brune PetscFunctionBegin; 16085fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1609422a814eSBarry Smith linesearch->result = result; 16106a388c36SPeter Brune PetscFunctionReturn(0); 16116a388c36SPeter Brune } 16126a388c36SPeter Brune 16139bd66eb0SPeter Brune /*@C 1614f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 16159bd66eb0SPeter Brune 16169bd66eb0SPeter Brune Input Parameters: 16179bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 16189bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 16199bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16209bd66eb0SPeter Brune 16219bd66eb0SPeter Brune Logically Collective on SNES 16229bd66eb0SPeter Brune 16239bd66eb0SPeter Brune Calling sequence of projectfunc: 16249bd66eb0SPeter Brune .vb 16259bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 16269bd66eb0SPeter Brune .ve 16279bd66eb0SPeter Brune 16289bd66eb0SPeter Brune Input parameters for projectfunc: 16299bd66eb0SPeter Brune + snes - nonlinear context 16309bd66eb0SPeter Brune - X - current solution 16319bd66eb0SPeter Brune 1632cd7522eaSPeter Brune Output parameters for projectfunc: 16339bd66eb0SPeter Brune . X - Projected solution 16349bd66eb0SPeter Brune 16359bd66eb0SPeter Brune Calling sequence of normfunc: 16369bd66eb0SPeter Brune .vb 16379bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 16389bd66eb0SPeter Brune .ve 16399bd66eb0SPeter Brune 1640cd7522eaSPeter Brune Input parameters for normfunc: 16419bd66eb0SPeter Brune + snes - nonlinear context 16429bd66eb0SPeter Brune . X - current solution 16439bd66eb0SPeter Brune - F - current residual 16449bd66eb0SPeter Brune 1645cd7522eaSPeter Brune Output parameters for normfunc: 16469bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 16479bd66eb0SPeter Brune 1648cd7522eaSPeter Brune Notes: 1649cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1650cd7522eaSPeter Brune 1651cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1652cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 16539bd66eb0SPeter Brune 16549bd66eb0SPeter Brune Level: developer 16559bd66eb0SPeter Brune 1656db781477SPatrick Sanan .seealso: `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 16579bd66eb0SPeter Brune @*/ 165825acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 16599bd66eb0SPeter Brune { 16609bd66eb0SPeter Brune PetscFunctionBegin; 1661f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 16629bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 16639bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 16649bd66eb0SPeter Brune PetscFunctionReturn(0); 16659bd66eb0SPeter Brune } 16669bd66eb0SPeter Brune 16679bd66eb0SPeter Brune /*@C 1668f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 16699bd66eb0SPeter Brune 1670f899ff85SJose E. Roman Input Parameter: 1671907376e6SBarry Smith . linesearch - the line search context, obtain with SNESGetLineSearch() 16729bd66eb0SPeter Brune 16739bd66eb0SPeter Brune Output Parameters: 16749bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16759bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16769bd66eb0SPeter Brune 16779bd66eb0SPeter Brune Logically Collective on SNES 16789bd66eb0SPeter Brune 16799bd66eb0SPeter Brune Level: developer 16809bd66eb0SPeter Brune 1681db781477SPatrick Sanan .seealso: `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 16829bd66eb0SPeter Brune @*/ 168325acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 16849bd66eb0SPeter Brune { 16859bd66eb0SPeter Brune PetscFunctionBegin; 16869bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16879bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16889bd66eb0SPeter Brune PetscFunctionReturn(0); 16899bd66eb0SPeter Brune } 16909bd66eb0SPeter Brune 1691bf7f4e0aSPeter Brune /*@C 16921c84c290SBarry Smith SNESLineSearchRegister - See SNESLineSearchRegister() 1693bf7f4e0aSPeter Brune 1694bf7f4e0aSPeter Brune Level: advanced 1695bf7f4e0aSPeter Brune @*/ 1696bdf89e91SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1697bf7f4e0aSPeter Brune { 1698bf7f4e0aSPeter Brune PetscFunctionBegin; 16999566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage()); 17009566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList,sname,function)); 1701bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1702bf7f4e0aSPeter Brune } 1703