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; 7f1c6b773SPeter Brune PetscLogEvent SNESLineSearch_Apply; 8bf7f4e0aSPeter Brune 9bf7f4e0aSPeter Brune #undef __FUNCT__ 10f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate" 11f40b411bSPeter Brune /*@ 12cd7522eaSPeter Brune SNESLineSearchCreate - Creates the line search context. 13f40b411bSPeter Brune 14cd7522eaSPeter Brune Logically Collective on Comm 15f40b411bSPeter Brune 16f40b411bSPeter Brune Input Parameters: 17cd7522eaSPeter Brune . comm - MPI communicator for the line search (typically from the associated SNES context). 18f40b411bSPeter Brune 19f40b411bSPeter Brune Output Parameters: 208e557f58SPeter Brune . outlinesearch - the new linesearch context 21f40b411bSPeter Brune 22162e0bf5SPeter Brune Level: developer 23162e0bf5SPeter Brune 24162e0bf5SPeter Brune Notes: 25162e0bf5SPeter Brune The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance 26162e0bf5SPeter Brune already associated with the SNES. This function is for developer use. 27f40b411bSPeter Brune 28cd7522eaSPeter Brune .keywords: LineSearch, create, context 29f40b411bSPeter Brune 30162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch() 31f40b411bSPeter Brune @*/ 32f40b411bSPeter Brune 33bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 34bf388a1fSBarry Smith { 35bf7f4e0aSPeter Brune PetscErrorCode ierr; 36f1c6b773SPeter Brune SNESLineSearch linesearch; 37bf388a1fSBarry Smith 38bf7f4e0aSPeter Brune PetscFunctionBegin; 39ea5d4fccSPeter Brune PetscValidPointer(outlinesearch,2); 40f34a81feSMatthew G. Knepley ierr = SNESInitializePackage();CHKERRQ(ierr); 410298fd71SBarry Smith *outlinesearch = NULL; 42f5af7f23SKarl Rupp 4373107ff1SLisandro Dalcin ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr); 44bf7f4e0aSPeter Brune 450298fd71SBarry Smith linesearch->vec_sol_new = NULL; 460298fd71SBarry Smith linesearch->vec_func_new = NULL; 470298fd71SBarry Smith linesearch->vec_sol = NULL; 480298fd71SBarry Smith linesearch->vec_func = NULL; 490298fd71SBarry Smith linesearch->vec_update = NULL; 509bd66eb0SPeter Brune 51bf7f4e0aSPeter Brune linesearch->lambda = 1.0; 52bf7f4e0aSPeter Brune linesearch->fnorm = 1.0; 53bf7f4e0aSPeter Brune linesearch->ynorm = 1.0; 54bf7f4e0aSPeter Brune linesearch->xnorm = 1.0; 55422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 56bf7f4e0aSPeter Brune linesearch->norms = PETSC_TRUE; 57bf7f4e0aSPeter Brune linesearch->keeplambda = PETSC_FALSE; 58bf7f4e0aSPeter Brune linesearch->damping = 1.0; 59bf7f4e0aSPeter Brune linesearch->maxstep = 1e8; 60bf7f4e0aSPeter Brune linesearch->steptol = 1e-12; 61516fe3c3SPeter Brune linesearch->rtol = 1e-8; 62516fe3c3SPeter Brune linesearch->atol = 1e-15; 63516fe3c3SPeter Brune linesearch->ltol = 1e-8; 640298fd71SBarry Smith linesearch->precheckctx = NULL; 650298fd71SBarry Smith linesearch->postcheckctx = NULL; 6659405d5eSPeter Brune linesearch->max_its = 1; 67bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 68bf7f4e0aSPeter Brune *outlinesearch = linesearch; 69bf7f4e0aSPeter Brune PetscFunctionReturn(0); 70bf7f4e0aSPeter Brune } 71bf7f4e0aSPeter Brune 72bf7f4e0aSPeter Brune #undef __FUNCT__ 73f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp" 74f40b411bSPeter Brune /*@ 7578bcb3b5SPeter Brune SNESLineSearchSetUp - Prepares the line search for being applied by allocating 7678bcb3b5SPeter Brune any required vectors. 77f40b411bSPeter Brune 78cd7522eaSPeter Brune Collective on SNESLineSearch 79f40b411bSPeter Brune 80f40b411bSPeter Brune Input Parameters: 81f40b411bSPeter Brune . linesearch - The LineSearch instance. 82f40b411bSPeter Brune 83cd7522eaSPeter Brune Notes: 84f190f2fcSBarry Smith For most cases, this needn't be called by users or outside of SNESLineSearchApply(). 85cd7522eaSPeter Brune The only current case where this is called outside of this is for the VI 8678bcb3b5SPeter Brune solvers, which modify the solution and work vectors before the first call 87cd7522eaSPeter Brune of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be 88cd7522eaSPeter Brune allocated upfront. 89cd7522eaSPeter Brune 9078bcb3b5SPeter Brune Level: advanced 91f40b411bSPeter Brune 92f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp 93f40b411bSPeter Brune 94f1c6b773SPeter Brune .seealso: SNESLineSearchReset() 95f40b411bSPeter Brune @*/ 96f40b411bSPeter Brune 97bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 98bf388a1fSBarry Smith { 99bf7f4e0aSPeter Brune PetscErrorCode ierr; 100bf388a1fSBarry Smith 101bf7f4e0aSPeter Brune PetscFunctionBegin; 102bf7f4e0aSPeter Brune if (!((PetscObject)linesearch)->type_name) { 1031a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 104bf7f4e0aSPeter Brune } 105bf7f4e0aSPeter Brune if (!linesearch->setupcalled) { 1069bd66eb0SPeter Brune if (!linesearch->vec_sol_new) { 107bf7f4e0aSPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr); 1089bd66eb0SPeter Brune } 1099bd66eb0SPeter Brune if (!linesearch->vec_func_new) { 1109bd66eb0SPeter Brune ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr); 1119bd66eb0SPeter Brune } 112bf7f4e0aSPeter Brune if (linesearch->ops->setup) { 113bf7f4e0aSPeter Brune ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr); 114bf7f4e0aSPeter Brune } 115ed07d7d7SPeter Brune if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);} 116bf7f4e0aSPeter Brune linesearch->lambda = linesearch->damping; 117bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_TRUE; 118bf7f4e0aSPeter Brune } 119bf7f4e0aSPeter Brune PetscFunctionReturn(0); 120bf7f4e0aSPeter Brune } 121bf7f4e0aSPeter Brune 122bf7f4e0aSPeter Brune #undef __FUNCT__ 123f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset" 124f40b411bSPeter Brune 125f40b411bSPeter Brune /*@ 126f190f2fcSBarry Smith SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search. 127f40b411bSPeter Brune 128f1c6b773SPeter Brune Collective on SNESLineSearch 129f40b411bSPeter Brune 130f40b411bSPeter Brune Input Parameters: 131f40b411bSPeter Brune . linesearch - The LineSearch instance. 132f40b411bSPeter Brune 133f190f2fcSBarry Smith Notes: Usually only called by SNESReset() 134f190f2fcSBarry Smith 135f190f2fcSBarry Smith Level: developer 136f40b411bSPeter Brune 137cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset 138f40b411bSPeter Brune 139f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp() 140f40b411bSPeter Brune @*/ 141f40b411bSPeter Brune 142bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 143bf388a1fSBarry Smith { 144bf7f4e0aSPeter Brune PetscErrorCode ierr; 145bf388a1fSBarry Smith 146bf7f4e0aSPeter Brune PetscFunctionBegin; 147f5af7f23SKarl Rupp if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch); 148f5af7f23SKarl Rupp 149bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr); 150bf7f4e0aSPeter Brune ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr); 151bf7f4e0aSPeter Brune 152bf7f4e0aSPeter Brune ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr); 153f5af7f23SKarl Rupp 154bf7f4e0aSPeter Brune linesearch->nwork = 0; 155bf7f4e0aSPeter Brune linesearch->setupcalled = PETSC_FALSE; 156bf7f4e0aSPeter Brune PetscFunctionReturn(0); 157bf7f4e0aSPeter Brune } 158bf7f4e0aSPeter Brune 159ed07d7d7SPeter Brune #undef __FUNCT__ 160ed07d7d7SPeter Brune #define __FUNCT__ "SNESLineSearchSetFunction" 161ed07d7d7SPeter Brune /*@C 162f190f2fcSBarry Smith SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search 163ed07d7d7SPeter Brune 164ed07d7d7SPeter Brune Input Parameters: 165ed07d7d7SPeter Brune . linesearch - the SNESLineSearch context 166f190f2fcSBarry Smith + func - function evaluation routine 167ed07d7d7SPeter Brune 168ed07d7d7SPeter Brune Level: developer 169ed07d7d7SPeter Brune 170f190f2fcSBarry Smith Notes: This is used internally by PETSc and not called by users 171f190f2fcSBarry Smith 172ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check 173ed07d7d7SPeter Brune 174f190f2fcSBarry Smith .seealso: SNESSetFunction() 175ed07d7d7SPeter Brune @*/ 176ed07d7d7SPeter Brune PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec)) 177ed07d7d7SPeter Brune { 178ed07d7d7SPeter Brune PetscFunctionBegin; 179ed07d7d7SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 180ed07d7d7SPeter Brune linesearch->ops->snesfunc = func; 181ed07d7d7SPeter Brune PetscFunctionReturn(0); 182ed07d7d7SPeter Brune } 183ed07d7d7SPeter Brune 184ed07d7d7SPeter Brune 1856b2b7091SBarry Smith /*MC 186f190f2fcSBarry Smith SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called 1876b2b7091SBarry Smith 1886b2b7091SBarry Smith Synopsis: 189aaa7dc30SBarry Smith #include <petscsnes.h> 1906b2b7091SBarry Smith SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed); 1916b2b7091SBarry Smith 1926b2b7091SBarry Smith Input Parameters: 1936b2b7091SBarry Smith + x - solution vector 1946b2b7091SBarry Smith . y - search direction vector 1956b2b7091SBarry Smith - changed - flag to indicate the precheck changed x or y. 1966b2b7091SBarry Smith 197f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck() 198f190f2fcSBarry Smith and SNESLineSearchGetPreCheck() 199f190f2fcSBarry Smith 200878cb397SSatish Balay Level: advanced 201878cb397SSatish Balay 202f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck() 2036b2b7091SBarry Smith M*/ 2046b2b7091SBarry Smith 20586d74e61SPeter Brune #undef __FUNCT__ 206f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck" 20786d74e61SPeter Brune /*@C 208f190f2fcSBarry Smith SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 209df3898eeSBarry Smith before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 210f190f2fcSBarry Smith determined the search direction. 21186d74e61SPeter Brune 212f1c6b773SPeter Brune Logically Collective on SNESLineSearch 21386d74e61SPeter Brune 21486d74e61SPeter Brune Input Parameters: 215f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 216f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence 217f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 21886d74e61SPeter Brune 21986d74e61SPeter Brune Level: intermediate 22086d74e61SPeter Brune 22186d74e61SPeter Brune .keywords: set, linesearch, pre-check 22286d74e61SPeter Brune 223f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 22486d74e61SPeter Brune @*/ 225f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx) 22686d74e61SPeter Brune { 2279bd66eb0SPeter Brune PetscFunctionBegin; 228f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 229f190f2fcSBarry Smith if (func) linesearch->ops->precheck = func; 23086d74e61SPeter Brune if (ctx) linesearch->precheckctx = ctx; 23186d74e61SPeter Brune PetscFunctionReturn(0); 23286d74e61SPeter Brune } 23386d74e61SPeter Brune 23486d74e61SPeter Brune #undef __FUNCT__ 235f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck" 23686d74e61SPeter Brune /*@C 23753deb899SBarry Smith SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 23886d74e61SPeter Brune 23986d74e61SPeter Brune Input Parameters: 240f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 24186d74e61SPeter Brune 24286d74e61SPeter Brune Output Parameters: 243f190f2fcSBarry Smith + func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence 244f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 24586d74e61SPeter Brune 24686d74e61SPeter Brune Level: intermediate 24786d74e61SPeter Brune 24886d74e61SPeter Brune .keywords: get, linesearch, pre-check 24986d74e61SPeter Brune 250f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck() 25186d74e61SPeter Brune @*/ 252f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx) 25386d74e61SPeter Brune { 2549bd66eb0SPeter Brune PetscFunctionBegin; 255f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 256f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck; 25786d74e61SPeter Brune if (ctx) *ctx = linesearch->precheckctx; 25886d74e61SPeter Brune PetscFunctionReturn(0); 25986d74e61SPeter Brune } 26086d74e61SPeter Brune 2616b2b7091SBarry Smith /*MC 262f190f2fcSBarry Smith SNESLineSearchPostCheckFunction - form of function that is called after line search is complete 2636b2b7091SBarry Smith 2646b2b7091SBarry Smith Synopsis: 265aaa7dc30SBarry Smith #include <petscsnes.h> 2666b2b7091SBarry Smith SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y, Vec w, *changed_y, PetscBool *changed_w); 2676b2b7091SBarry Smith 2686b2b7091SBarry Smith Input Parameters: 2696b2b7091SBarry Smith + x - old solution vector 2706b2b7091SBarry Smith . y - search direction vector 2716b2b7091SBarry Smith . w - new solution vector 2726b2b7091SBarry Smith . changed_y - indicates that the line search changed y 2736b2b7091SBarry Smith - changed_w - indicates that the line search changed w 2746b2b7091SBarry Smith 275f190f2fcSBarry Smith Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck() 276f190f2fcSBarry Smith and SNESLineSearchGetPostCheck() 277f190f2fcSBarry Smith 278878cb397SSatish Balay Level: advanced 2796b2b7091SBarry Smith 280f190f2fcSBarry Smith .seealso: SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck() 2816b2b7091SBarry Smith M*/ 2826b2b7091SBarry Smith 28386d74e61SPeter Brune #undef __FUNCT__ 284f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck" 28586d74e61SPeter Brune /*@C 286f190f2fcSBarry Smith SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 287f190f2fcSBarry Smith direction and length. Allows the user a chance to change or override the decision of the line search routine 28886d74e61SPeter Brune 289f1c6b773SPeter Brune Logically Collective on SNESLineSearch 29086d74e61SPeter Brune 29186d74e61SPeter Brune Input Parameters: 292f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 293f190f2fcSBarry Smith . func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence 294f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 29586d74e61SPeter Brune 29686d74e61SPeter Brune Level: intermediate 29786d74e61SPeter Brune 29886d74e61SPeter Brune .keywords: set, linesearch, post-check 29986d74e61SPeter Brune 300f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck() 30186d74e61SPeter Brune @*/ 302f190f2fcSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 30386d74e61SPeter Brune { 30486d74e61SPeter Brune PetscFunctionBegin; 305f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 306f190f2fcSBarry Smith if (func) linesearch->ops->postcheck = func; 30786d74e61SPeter Brune if (ctx) linesearch->postcheckctx = ctx; 30886d74e61SPeter Brune PetscFunctionReturn(0); 30986d74e61SPeter Brune } 31086d74e61SPeter Brune 31186d74e61SPeter Brune #undef __FUNCT__ 312f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck" 31386d74e61SPeter Brune /*@C 314f1c6b773SPeter Brune SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 31586d74e61SPeter Brune 31686d74e61SPeter Brune Input Parameters: 317f1c6b773SPeter Brune . linesearch - the SNESLineSearch context 31886d74e61SPeter Brune 31986d74e61SPeter Brune Output Parameters: 320f190f2fcSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction 321f190f2fcSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 32286d74e61SPeter Brune 32386d74e61SPeter Brune Level: intermediate 32486d74e61SPeter Brune 32586d74e61SPeter Brune .keywords: get, linesearch, post-check 32686d74e61SPeter Brune 327f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck() 32886d74e61SPeter Brune @*/ 329f190f2fcSBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 33086d74e61SPeter Brune { 3319bd66eb0SPeter Brune PetscFunctionBegin; 332f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 333f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck; 33486d74e61SPeter Brune if (ctx) *ctx = linesearch->postcheckctx; 33586d74e61SPeter Brune PetscFunctionReturn(0); 33686d74e61SPeter Brune } 33786d74e61SPeter Brune 338bf7f4e0aSPeter Brune #undef __FUNCT__ 339f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck" 340f40b411bSPeter Brune /*@ 341f1c6b773SPeter Brune SNESLineSearchPreCheck - Prepares the line search for being applied. 342f40b411bSPeter Brune 343cd7522eaSPeter Brune Logically Collective on SNESLineSearch 344f40b411bSPeter Brune 345f40b411bSPeter Brune Input Parameters: 3467b1df9c1SPeter Brune + linesearch - The linesearch instance. 3477b1df9c1SPeter Brune . X - The current solution 3487b1df9c1SPeter Brune - Y - The step direction 349f40b411bSPeter Brune 350f40b411bSPeter Brune Output Parameters: 3518e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything 352f40b411bSPeter Brune 353f190f2fcSBarry Smith Level: developer 354f40b411bSPeter Brune 355f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 356f40b411bSPeter Brune 357f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck() 358f40b411bSPeter Brune @*/ 3597b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 360bf7f4e0aSPeter Brune { 361bf7f4e0aSPeter Brune PetscErrorCode ierr; 3625fd66863SKarl Rupp 363bf7f4e0aSPeter Brune PetscFunctionBegin; 364bf7f4e0aSPeter Brune *changed = PETSC_FALSE; 3656b2b7091SBarry Smith if (linesearch->ops->precheck) { 3666b2b7091SBarry Smith ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr); 36738bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed,4); 368bf7f4e0aSPeter Brune } 369bf7f4e0aSPeter Brune PetscFunctionReturn(0); 370bf7f4e0aSPeter Brune } 371bf7f4e0aSPeter Brune 372bf7f4e0aSPeter Brune #undef __FUNCT__ 373f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck" 374f40b411bSPeter Brune /*@ 375f1c6b773SPeter Brune SNESLineSearchPostCheck - Prepares the line search for being applied. 376f40b411bSPeter Brune 377cd7522eaSPeter Brune Logically Collective on SNESLineSearch 378f40b411bSPeter Brune 379f40b411bSPeter Brune Input Parameters: 3807b1df9c1SPeter Brune + linesearch - The linesearch context 3817b1df9c1SPeter Brune . X - The last solution 3827b1df9c1SPeter Brune . Y - The step direction 3837b1df9c1SPeter Brune - W - The updated solution, W = X + lambda*Y for some lambda 384f40b411bSPeter Brune 385f40b411bSPeter Brune Output Parameters: 38678bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed. 38778bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed. 388f40b411bSPeter Brune 389f190f2fcSBarry Smith Level: developer 390f40b411bSPeter Brune 391f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 392f40b411bSPeter Brune 393f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck() 394f40b411bSPeter Brune @*/ 3957b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 396bf7f4e0aSPeter Brune { 397bf7f4e0aSPeter Brune PetscErrorCode ierr; 398bf388a1fSBarry Smith 399bf7f4e0aSPeter Brune PetscFunctionBegin; 400bf7f4e0aSPeter Brune *changed_Y = PETSC_FALSE; 401bf7f4e0aSPeter Brune *changed_W = PETSC_FALSE; 4026b2b7091SBarry Smith if (linesearch->ops->postcheck) { 4036b2b7091SBarry Smith ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr); 40438bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5); 40538bcdd5aSPeter Brune PetscValidLogicalCollectiveBool(linesearch,*changed_W,6); 40686d74e61SPeter Brune } 40786d74e61SPeter Brune PetscFunctionReturn(0); 40886d74e61SPeter Brune } 40986d74e61SPeter Brune 41086d74e61SPeter Brune #undef __FUNCT__ 411f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard" 41286d74e61SPeter Brune /*@C 41386d74e61SPeter Brune SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 41486d74e61SPeter Brune 415cd7522eaSPeter Brune Logically Collective on SNESLineSearch 41686d74e61SPeter Brune 41786d74e61SPeter Brune Input Arguments: 41886d74e61SPeter Brune + linesearch - linesearch context 41986d74e61SPeter Brune . X - base state for this step 42086d74e61SPeter Brune . Y - initial correction 421907376e6SBarry Smith - ctx - context for this function 42286d74e61SPeter Brune 42386d74e61SPeter Brune Output Arguments: 42486d74e61SPeter Brune + Y - correction, possibly modified 42586d74e61SPeter Brune - changed - flag indicating that Y was modified 42686d74e61SPeter Brune 42786d74e61SPeter Brune Options Database Key: 428cd7522eaSPeter Brune + -snes_linesearch_precheck_picard - activate this routine 429cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle 43086d74e61SPeter Brune 43186d74e61SPeter Brune Level: advanced 43286d74e61SPeter Brune 43386d74e61SPeter Brune Notes: 43486d74e61SPeter Brune This function should be passed to SNESLineSearchSetPreCheck() 43586d74e61SPeter Brune 43686d74e61SPeter Brune The justification for this method involves the linear convergence of a Picard iteration 43786d74e61SPeter Brune so the Picard linearization should be provided in place of the "Jacobian". This correction 43886d74e61SPeter Brune is generally not useful when using a Newton linearization. 43986d74e61SPeter Brune 44086d74e61SPeter Brune Reference: 44186d74e61SPeter Brune Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 44286d74e61SPeter Brune 44386d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck() 44486d74e61SPeter Brune @*/ 445f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx) 44686d74e61SPeter Brune { 44786d74e61SPeter Brune PetscErrorCode ierr; 44886d74e61SPeter Brune PetscReal angle = *(PetscReal*)linesearch->precheckctx; 44986d74e61SPeter Brune Vec Ylast; 45086d74e61SPeter Brune PetscScalar dot; 45186d74e61SPeter Brune PetscInt iter; 45286d74e61SPeter Brune PetscReal ynorm,ylastnorm,theta,angle_radians; 45386d74e61SPeter Brune SNES snes; 45486d74e61SPeter Brune 45586d74e61SPeter Brune PetscFunctionBegin; 456f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 45786d74e61SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 45886d74e61SPeter Brune if (!Ylast) { 45986d74e61SPeter Brune ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 46086d74e61SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 46186d74e61SPeter Brune ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 46286d74e61SPeter Brune } 46386d74e61SPeter Brune ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 46486d74e61SPeter Brune if (iter < 2) { 46586d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 46686d74e61SPeter Brune *changed = PETSC_FALSE; 46786d74e61SPeter Brune PetscFunctionReturn(0); 46886d74e61SPeter Brune } 46986d74e61SPeter Brune 47086d74e61SPeter Brune ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 47186d74e61SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 47286d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 47386d74e61SPeter Brune /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 474255453a1SBarry Smith theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 47586d74e61SPeter Brune angle_radians = angle * PETSC_PI / 180.; 47686d74e61SPeter Brune if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 47786d74e61SPeter Brune /* Modify the step Y */ 47886d74e61SPeter Brune PetscReal alpha,ydiffnorm; 47986d74e61SPeter Brune ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 48086d74e61SPeter Brune ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 48186d74e61SPeter Brune alpha = ylastnorm / ydiffnorm; 48286d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 48386d74e61SPeter Brune ierr = VecScale(Y,alpha);CHKERRQ(ierr); 484c69d1a72SBarry Smith ierr = PetscInfo3(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);CHKERRQ(ierr); 48586d74e61SPeter Brune } else { 486c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr); 48786d74e61SPeter Brune ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 48886d74e61SPeter Brune *changed = PETSC_FALSE; 489bf7f4e0aSPeter Brune } 490bf7f4e0aSPeter Brune PetscFunctionReturn(0); 491bf7f4e0aSPeter Brune } 492bf7f4e0aSPeter Brune 493bf7f4e0aSPeter Brune #undef __FUNCT__ 494f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply" 495f40b411bSPeter Brune /*@ 496cd7522eaSPeter Brune SNESLineSearchApply - Computes the line-search update. 497f40b411bSPeter Brune 498f1c6b773SPeter Brune Collective on SNESLineSearch 499f40b411bSPeter Brune 500f40b411bSPeter Brune Input Parameters: 5018e557f58SPeter Brune + linesearch - The linesearch context 5028e557f58SPeter Brune . X - The current solution 5038e557f58SPeter Brune . F - The current function 5048e557f58SPeter Brune . fnorm - The current norm 5058e557f58SPeter Brune - Y - The search direction 506f40b411bSPeter Brune 507f40b411bSPeter Brune Output Parameters: 5088e557f58SPeter Brune + X - The new solution 5098e557f58SPeter Brune . F - The new function 5108e557f58SPeter Brune - fnorm - The new function norm 511f40b411bSPeter Brune 512cd7522eaSPeter Brune Options Database Keys: 513d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell 514cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 5153c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 516cd7522eaSPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms 5173c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 5183c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches 519cd7522eaSPeter Brune 520cd7522eaSPeter Brune Notes: 521cd7522eaSPeter Brune This is typically called from within a SNESSolve() implementation in order to 522cd7522eaSPeter Brune help with convergence of the nonlinear method. Various SNES types use line searches 523cd7522eaSPeter Brune in different ways, but the overarching theme is that a line search is used to determine 524cd7522eaSPeter Brune an optimal damping parameter of a step at each iteration of the method. Each 525cd7522eaSPeter Brune application of the line search may invoke SNESComputeFunction several times, and 526cd7522eaSPeter Brune therefore may be fairly expensive. 527cd7522eaSPeter Brune 528f40b411bSPeter Brune Level: Intermediate 529f40b411bSPeter Brune 530f1c6b773SPeter Brune .keywords: SNESLineSearch, Create 531f40b411bSPeter Brune 532cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction() 533f40b411bSPeter Brune @*/ 534bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) 535bf388a1fSBarry Smith { 536bf7f4e0aSPeter Brune PetscErrorCode ierr; 537bf7f4e0aSPeter Brune 538bf388a1fSBarry Smith PetscFunctionBegin; 539f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 540bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 541bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 542bf7f4e0aSPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 543bf7f4e0aSPeter Brune 544422a814eSBarry Smith linesearch->result = SNES_LINESEARCH_SUCCEEDED; 545bf7f4e0aSPeter Brune 546bf7f4e0aSPeter Brune linesearch->vec_sol = X; 547bf7f4e0aSPeter Brune linesearch->vec_update = Y; 548bf7f4e0aSPeter Brune linesearch->vec_func = F; 549bf7f4e0aSPeter Brune 550f1c6b773SPeter Brune ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr); 551bf7f4e0aSPeter Brune 552f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 553bf7f4e0aSPeter Brune 554f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm; 555f5af7f23SKarl Rupp else { 556bf7f4e0aSPeter Brune ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 557bf7f4e0aSPeter Brune } 558bf7f4e0aSPeter Brune 559f1c6b773SPeter Brune ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 560bf7f4e0aSPeter Brune 561bf7f4e0aSPeter Brune ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 562bf7f4e0aSPeter Brune 563f1c6b773SPeter Brune ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 564bf7f4e0aSPeter Brune 565f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 566bf7f4e0aSPeter Brune PetscFunctionReturn(0); 567bf7f4e0aSPeter Brune } 568bf7f4e0aSPeter Brune 569bf7f4e0aSPeter Brune #undef __FUNCT__ 570f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy" 571f40b411bSPeter Brune /*@ 572f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance. 573f40b411bSPeter Brune 574f1c6b773SPeter Brune Collective on SNESLineSearch 575f40b411bSPeter Brune 576f40b411bSPeter Brune Input Parameters: 5778e557f58SPeter Brune . linesearch - The linesearch context 578f40b411bSPeter Brune 579f40b411bSPeter Brune Level: Intermediate 580f40b411bSPeter Brune 58178bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy 582f40b411bSPeter Brune 583cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy() 584f40b411bSPeter Brune @*/ 585bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) 586bf388a1fSBarry Smith { 587bf7f4e0aSPeter Brune PetscErrorCode ierr; 588bf388a1fSBarry Smith 589bf7f4e0aSPeter Brune PetscFunctionBegin; 590bf7f4e0aSPeter Brune if (!*linesearch) PetscFunctionReturn(0); 591f1c6b773SPeter Brune PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1); 592bf7f4e0aSPeter Brune if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 593e04113cfSBarry Smith ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr); 59422d28d08SBarry Smith ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr); 595f5af7f23SKarl Rupp if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch); 596bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 597e7058c64SPeter Brune ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 598bf7f4e0aSPeter Brune PetscFunctionReturn(0); 599bf7f4e0aSPeter Brune } 600bf7f4e0aSPeter Brune 601bf7f4e0aSPeter Brune #undef __FUNCT__ 602f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor" 603f40b411bSPeter Brune /*@ 604cd7522eaSPeter Brune SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search. 605bf7f4e0aSPeter Brune 606bf7f4e0aSPeter Brune Input Parameters: 607bf7f4e0aSPeter Brune + snes - nonlinear context obtained from SNESCreate() 608bf7f4e0aSPeter Brune - flg - PETSC_TRUE to monitor the line search 609bf7f4e0aSPeter Brune 610bf7f4e0aSPeter Brune Logically Collective on SNES 611bf7f4e0aSPeter Brune 612bf7f4e0aSPeter Brune Options Database: 6138e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 614bf7f4e0aSPeter Brune 615bf7f4e0aSPeter Brune Level: intermediate 616bf7f4e0aSPeter Brune 617cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer 618bf7f4e0aSPeter Brune @*/ 619f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg) 620bf7f4e0aSPeter Brune { 621bf7f4e0aSPeter Brune PetscErrorCode ierr; 622bf388a1fSBarry Smith 623bf7f4e0aSPeter Brune PetscFunctionBegin; 624bf7f4e0aSPeter Brune if (flg && !linesearch->monitor) { 625ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr); 626bf7f4e0aSPeter Brune } else if (!flg && linesearch->monitor) { 627bf7f4e0aSPeter Brune ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 628bf7f4e0aSPeter Brune } 629bf7f4e0aSPeter Brune PetscFunctionReturn(0); 630bf7f4e0aSPeter Brune } 631bf7f4e0aSPeter Brune 632bf7f4e0aSPeter Brune #undef __FUNCT__ 633f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor" 634f40b411bSPeter Brune /*@ 635cd7522eaSPeter Brune SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor. 6366a388c36SPeter Brune 637f190f2fcSBarry Smith Input Parameter: 6388e557f58SPeter Brune . linesearch - linesearch context 639f40b411bSPeter Brune 640f190f2fcSBarry Smith Output Parameter: 6418e557f58SPeter Brune . monitor - monitor context 642f40b411bSPeter Brune 643f40b411bSPeter Brune Logically Collective on SNES 644f40b411bSPeter Brune 645f40b411bSPeter Brune Options Database Keys: 6468e557f58SPeter Brune . -snes_linesearch_monitor - enables the monitor 647f40b411bSPeter Brune 648f40b411bSPeter Brune Level: intermediate 649f40b411bSPeter Brune 650cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer 651f40b411bSPeter Brune @*/ 652f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 6536a388c36SPeter Brune { 6546a388c36SPeter Brune PetscFunctionBegin; 655f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 6566a388c36SPeter Brune if (monitor) { 6576a388c36SPeter Brune PetscValidPointer(monitor, 2); 6586a388c36SPeter Brune *monitor = linesearch->monitor; 6596a388c36SPeter Brune } 6606a388c36SPeter Brune PetscFunctionReturn(0); 6616a388c36SPeter Brune } 6626a388c36SPeter Brune 6636a388c36SPeter Brune #undef __FUNCT__ 664f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions" 665f40b411bSPeter Brune /*@ 666f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search 667f40b411bSPeter Brune 668f40b411bSPeter Brune Input Parameters: 6698e557f58SPeter Brune . linesearch - linesearch context 670f40b411bSPeter Brune 671f40b411bSPeter Brune Options Database Keys: 672d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell 6733c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 6745a9b6599SPeter Brune . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch type 67571eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length 6761a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size 6771a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches 6781a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches 6791a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 6801a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches 681cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches 6828e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter 683cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 6841a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 6851a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method 686f40b411bSPeter Brune 687f1c6b773SPeter Brune Logically Collective on SNESLineSearch 688f40b411bSPeter Brune 689f40b411bSPeter Brune Level: intermediate 690f40b411bSPeter Brune 6913c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard() 692f40b411bSPeter Brune @*/ 693bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 694bf388a1fSBarry Smith { 695bf7f4e0aSPeter Brune PetscErrorCode ierr; 6961a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC; 697bf7f4e0aSPeter Brune char type[256]; 698bf7f4e0aSPeter Brune PetscBool flg, set; 699bf388a1fSBarry Smith 700bf7f4e0aSPeter Brune PetscFunctionBegin; 7010f51fdf8SToby Isaac ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr); 702bf7f4e0aSPeter Brune 703bf7f4e0aSPeter Brune ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 704f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 705a264d7a6SBarry Smith ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 706bf7f4e0aSPeter Brune if (flg) { 707f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 708bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) { 709f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 710bf7f4e0aSPeter Brune } 711bf7f4e0aSPeter Brune 7128afaa268SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 713f1c6b773SPeter Brune if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 714bf7f4e0aSPeter Brune 7151a9b3a06SPeter Brune /* tolerances */ 71694ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr); 71794ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr); 71894ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr); 71994ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr); 72094ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr); 72194ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr); 7227a35526eSPeter Brune 7231a9b3a06SPeter Brune /* damping parameters */ 72494ae4db5SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr); 7251a9b3a06SPeter Brune 72694ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr); 7271a9b3a06SPeter Brune 7281a9b3a06SPeter Brune /* precheck */ 7297a35526eSPeter Brune ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 7307a35526eSPeter Brune if (set) { 7317a35526eSPeter Brune if (flg) { 7327a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 733f5af7f23SKarl Rupp 7347a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 7350298fd71SBarry Smith "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr); 7367a35526eSPeter Brune ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 7377a35526eSPeter Brune } else { 7380298fd71SBarry Smith ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr); 7397a35526eSPeter Brune } 7407a35526eSPeter Brune } 74194ae4db5SBarry Smith ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr); 74294ae4db5SBarry Smith ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr); 7437a35526eSPeter Brune 7445a9b6599SPeter Brune if (linesearch->ops->setfromoptions) { 745e55864a3SBarry Smith (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr); 7465a9b6599SPeter Brune } 7475a9b6599SPeter Brune 748*0633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr); 749bf7f4e0aSPeter Brune ierr = PetscOptionsEnd();CHKERRQ(ierr); 750bf7f4e0aSPeter Brune PetscFunctionReturn(0); 751bf7f4e0aSPeter Brune } 752bf7f4e0aSPeter Brune 753bf7f4e0aSPeter Brune #undef __FUNCT__ 754f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView" 755f40b411bSPeter Brune /*@ 756f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search 757f40b411bSPeter Brune 758f40b411bSPeter Brune Input Parameters: 7598e557f58SPeter Brune . linesearch - linesearch context 760f40b411bSPeter Brune 761f1c6b773SPeter Brune Logically Collective on SNESLineSearch 762f40b411bSPeter Brune 763f40b411bSPeter Brune Level: intermediate 764f40b411bSPeter Brune 765f190f2fcSBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchMonitor() 766f40b411bSPeter Brune @*/ 767bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 768bf388a1fSBarry Smith { 7697f1410a3SPeter Brune PetscErrorCode ierr; 7707f1410a3SPeter Brune PetscBool iascii; 771bf388a1fSBarry Smith 772bf7f4e0aSPeter Brune PetscFunctionBegin; 7737f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 7747f1410a3SPeter Brune if (!viewer) { 775ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr); 7767f1410a3SPeter Brune } 7777f1410a3SPeter Brune PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 7787f1410a3SPeter Brune PetscCheckSameComm(linesearch,1,viewer,2); 779f40b411bSPeter Brune 780251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7817f1410a3SPeter Brune if (iascii) { 782dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr); 7837f1410a3SPeter Brune if (linesearch->ops->view) { 7847f1410a3SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7857f1410a3SPeter Brune ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 7867f1410a3SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7877f1410a3SPeter Brune } 788c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr); 789c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr); 7907f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 7916b2b7091SBarry Smith if (linesearch->ops->precheck) { 7926b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 7937f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 7947f1410a3SPeter Brune } else { 7957f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 7967f1410a3SPeter Brune } 7977f1410a3SPeter Brune } 7986b2b7091SBarry Smith if (linesearch->ops->postcheck) { 7997f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 8007f1410a3SPeter Brune } 8017f1410a3SPeter Brune } 802bf7f4e0aSPeter Brune PetscFunctionReturn(0); 803bf7f4e0aSPeter Brune } 804bf7f4e0aSPeter Brune 805bf7f4e0aSPeter Brune #undef __FUNCT__ 806f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType" 807ea5d4fccSPeter Brune /*@C 808f1c6b773SPeter Brune SNESLineSearchSetType - Sets the linesearch type 809f40b411bSPeter Brune 810f190f2fcSBarry Smith Logically Collective on SNESLineSearch 811f190f2fcSBarry Smith 812f40b411bSPeter Brune Input Parameters: 8138e557f58SPeter Brune + linesearch - linesearch context 814f40b411bSPeter Brune - type - The type of line search to be used 815f40b411bSPeter Brune 816cd7522eaSPeter Brune Available Types: 817cd7522eaSPeter Brune + basic - Simple damping line search. 8188e557f58SPeter Brune . bt - Backtracking line search over the L2 norm of the function 8198e557f58SPeter Brune . l2 - Secant line search over the L2 norm of the function 8208e557f58SPeter Brune . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 821d4c6564cSPatrick Farrell . nleqerr - Affine-covariant error-oriented linesearch 8228e557f58SPeter Brune - shell - User provided SNESLineSearch implementation 823cd7522eaSPeter Brune 824f40b411bSPeter Brune Level: intermediate 825f40b411bSPeter Brune 826f1c6b773SPeter Brune .seealso: SNESLineSearchCreate() 827f40b411bSPeter Brune @*/ 82819fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 829bf7f4e0aSPeter Brune { 830f1c6b773SPeter Brune PetscErrorCode ierr,(*r)(SNESLineSearch); 831bf7f4e0aSPeter Brune PetscBool match; 832bf7f4e0aSPeter Brune 833bf7f4e0aSPeter Brune PetscFunctionBegin; 834f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 835bf7f4e0aSPeter Brune PetscValidCharPointer(type,2); 836bf7f4e0aSPeter Brune 837251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 838bf7f4e0aSPeter Brune if (match) PetscFunctionReturn(0); 839bf7f4e0aSPeter Brune 8401c9cd337SJed Brown ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr); 841bf7f4e0aSPeter Brune if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 842bf7f4e0aSPeter Brune /* Destroy the previous private linesearch context */ 843bf7f4e0aSPeter Brune if (linesearch->ops->destroy) { 844bf7f4e0aSPeter Brune ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 845f5af7f23SKarl Rupp 8460298fd71SBarry Smith linesearch->ops->destroy = NULL; 847bf7f4e0aSPeter Brune } 848f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */ 849bf7f4e0aSPeter Brune linesearch->ops->apply = 0; 850bf7f4e0aSPeter Brune linesearch->ops->view = 0; 851bf7f4e0aSPeter Brune linesearch->ops->setfromoptions = 0; 852bf7f4e0aSPeter Brune linesearch->ops->destroy = 0; 853bf7f4e0aSPeter Brune 854bf7f4e0aSPeter Brune ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 855bf7f4e0aSPeter Brune ierr = (*r)(linesearch);CHKERRQ(ierr); 856bf7f4e0aSPeter Brune PetscFunctionReturn(0); 857bf7f4e0aSPeter Brune } 858bf7f4e0aSPeter Brune 859bf7f4e0aSPeter Brune #undef __FUNCT__ 860f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES" 861f40b411bSPeter Brune /*@ 86278bcb3b5SPeter Brune SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 863f40b411bSPeter Brune 864f40b411bSPeter Brune Input Parameters: 8658e557f58SPeter Brune + linesearch - linesearch context 866f40b411bSPeter Brune - snes - The snes instance 867f40b411bSPeter Brune 86878bcb3b5SPeter Brune Level: developer 86978bcb3b5SPeter Brune 87078bcb3b5SPeter Brune Notes: 871f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with 8727601faf0SJed Brown SNESGetLineSearch(). This routine is therefore mainly called within SNES 87378bcb3b5SPeter Brune implementations. 874f40b411bSPeter Brune 8758141a3b9SPeter Brune Level: developer 876f40b411bSPeter Brune 877cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 878f40b411bSPeter Brune @*/ 879bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 880bf388a1fSBarry Smith { 881bf7f4e0aSPeter Brune PetscFunctionBegin; 882f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 883bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 884bf7f4e0aSPeter Brune linesearch->snes = snes; 885bf7f4e0aSPeter Brune PetscFunctionReturn(0); 886bf7f4e0aSPeter Brune } 887bf7f4e0aSPeter Brune 888bf7f4e0aSPeter Brune #undef __FUNCT__ 889f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES" 890f40b411bSPeter Brune /*@ 8918141a3b9SPeter Brune SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 8928141a3b9SPeter Brune Having an associated SNES is necessary because most line search implementations must be able to 8938141a3b9SPeter Brune evaluate the function using SNESComputeFunction() for the associated SNES. This routine 8948141a3b9SPeter Brune is used in the line search implementations when one must get this associated SNES instance. 895f40b411bSPeter Brune 896f40b411bSPeter Brune Input Parameters: 8978e557f58SPeter Brune . linesearch - linesearch context 898f40b411bSPeter Brune 899f40b411bSPeter Brune Output Parameters: 900f40b411bSPeter Brune . snes - The snes instance 901f40b411bSPeter Brune 9028141a3b9SPeter Brune Level: developer 903f40b411bSPeter Brune 904cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 905f40b411bSPeter Brune @*/ 906bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 907bf388a1fSBarry Smith { 908bf7f4e0aSPeter Brune PetscFunctionBegin; 909f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9106a388c36SPeter Brune PetscValidPointer(snes, 2); 911bf7f4e0aSPeter Brune *snes = linesearch->snes; 912bf7f4e0aSPeter Brune PetscFunctionReturn(0); 913bf7f4e0aSPeter Brune } 914bf7f4e0aSPeter Brune 9156a388c36SPeter Brune #undef __FUNCT__ 916f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda" 917f40b411bSPeter Brune /*@ 918f1c6b773SPeter Brune SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 919f40b411bSPeter Brune 920f40b411bSPeter Brune Input Parameters: 9218e557f58SPeter Brune . linesearch - linesearch context 922f40b411bSPeter Brune 923f40b411bSPeter Brune Output Parameters: 924cd7522eaSPeter Brune . lambda - The last steplength computed during SNESLineSearchApply() 925f40b411bSPeter Brune 92678bcb3b5SPeter Brune Level: advanced 92778bcb3b5SPeter Brune 9288e557f58SPeter Brune Notes: 9298e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and 93078bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the 93178bcb3b5SPeter Brune solution and the function. For instance, SNESQN may be scaled by the 93278bcb3b5SPeter Brune line search lambda using the argument -snes_qn_scaling ls. 93378bcb3b5SPeter Brune 934cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 935f40b411bSPeter Brune @*/ 936f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 9376a388c36SPeter Brune { 9386a388c36SPeter Brune PetscFunctionBegin; 939f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9406a388c36SPeter Brune PetscValidPointer(lambda, 2); 9416a388c36SPeter Brune *lambda = linesearch->lambda; 9426a388c36SPeter Brune PetscFunctionReturn(0); 9436a388c36SPeter Brune } 9446a388c36SPeter Brune 9456a388c36SPeter Brune #undef __FUNCT__ 946f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda" 947f40b411bSPeter Brune /*@ 948f1c6b773SPeter Brune SNESLineSearchSetLambda - Sets the linesearch steplength. 949f40b411bSPeter Brune 950f40b411bSPeter Brune Input Parameters: 9518e557f58SPeter Brune + linesearch - linesearch context 952f40b411bSPeter Brune - lambda - The last steplength. 953f40b411bSPeter Brune 954cd7522eaSPeter Brune Notes: 955f190f2fcSBarry Smith This routine is typically used within implementations of SNESLineSearchApply() 956cd7522eaSPeter Brune to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 957cd7522eaSPeter Brune added in order to facilitate Quasi-Newton methods that use the previous steplength 958cd7522eaSPeter Brune as an inner scaling parameter. 959cd7522eaSPeter Brune 96078bcb3b5SPeter Brune Level: advanced 961f40b411bSPeter Brune 962f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda() 963f40b411bSPeter Brune @*/ 964f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 9656a388c36SPeter Brune { 9666a388c36SPeter Brune PetscFunctionBegin; 967f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 9686a388c36SPeter Brune linesearch->lambda = lambda; 9696a388c36SPeter Brune PetscFunctionReturn(0); 9706a388c36SPeter Brune } 9716a388c36SPeter Brune 9726a388c36SPeter Brune #undef __FUNCT__ 973f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances" 974f40b411bSPeter Brune /*@ 9753c7d6663SPeter Brune SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 97678bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 97778bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 97878bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 979f40b411bSPeter Brune 980f40b411bSPeter Brune Input Parameters: 9818e557f58SPeter Brune . linesearch - linesearch context 982f40b411bSPeter Brune 983f40b411bSPeter Brune Output Parameters: 984516fe3c3SPeter Brune + steptol - The minimum steplength 9856cc8e53bSPeter Brune . maxstep - The maximum steplength 986516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 987516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 988516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 989516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 990f40b411bSPeter Brune 99178bcb3b5SPeter Brune Level: intermediate 99278bcb3b5SPeter Brune 99378bcb3b5SPeter Brune Notes: 99478bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as 9953c7d6663SPeter Brune the type requires. 996516fe3c3SPeter Brune 997f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances() 998f40b411bSPeter Brune @*/ 999f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 10006a388c36SPeter Brune { 10016a388c36SPeter Brune PetscFunctionBegin; 1002f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1003516fe3c3SPeter Brune if (steptol) { 10046a388c36SPeter Brune PetscValidPointer(steptol, 2); 10056a388c36SPeter Brune *steptol = linesearch->steptol; 1006516fe3c3SPeter Brune } 1007516fe3c3SPeter Brune if (maxstep) { 1008516fe3c3SPeter Brune PetscValidPointer(maxstep, 3); 1009516fe3c3SPeter Brune *maxstep = linesearch->maxstep; 1010516fe3c3SPeter Brune } 1011516fe3c3SPeter Brune if (rtol) { 1012516fe3c3SPeter Brune PetscValidPointer(rtol, 4); 1013516fe3c3SPeter Brune *rtol = linesearch->rtol; 1014516fe3c3SPeter Brune } 1015516fe3c3SPeter Brune if (atol) { 1016516fe3c3SPeter Brune PetscValidPointer(atol, 5); 1017516fe3c3SPeter Brune *atol = linesearch->atol; 1018516fe3c3SPeter Brune } 1019516fe3c3SPeter Brune if (ltol) { 1020516fe3c3SPeter Brune PetscValidPointer(ltol, 6); 1021516fe3c3SPeter Brune *ltol = linesearch->ltol; 1022516fe3c3SPeter Brune } 1023516fe3c3SPeter Brune if (max_its) { 1024516fe3c3SPeter Brune PetscValidPointer(max_its, 7); 1025516fe3c3SPeter Brune *max_its = linesearch->max_its; 1026516fe3c3SPeter Brune } 10276a388c36SPeter Brune PetscFunctionReturn(0); 10286a388c36SPeter Brune } 10296a388c36SPeter Brune 10306a388c36SPeter Brune #undef __FUNCT__ 1031f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances" 1032f40b411bSPeter Brune /*@ 10333c7d6663SPeter Brune SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 103478bcb3b5SPeter Brune tolerances for the relative and absolute change in the function norm, the change 103578bcb3b5SPeter Brune in lambda for iterative line searches, the minimum steplength, the maximum steplength, 103678bcb3b5SPeter Brune and the maximum number of iterations the line search procedure may take. 1037f40b411bSPeter Brune 1038f40b411bSPeter Brune Input Parameters: 10398e557f58SPeter Brune + linesearch - linesearch context 1040516fe3c3SPeter Brune . steptol - The minimum steplength 10416cc8e53bSPeter Brune . maxstep - The maximum steplength 1042516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches 1043516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches 1044516fe3c3SPeter Brune . ltol - The change in lambda tolerance for iterative line searches 1045516fe3c3SPeter Brune - max_it - The maximum number of iterations of the line search 1046f40b411bSPeter Brune 104778bcb3b5SPeter Brune Notes: 10483c7d6663SPeter Brune The user may choose to not set any of the tolerances using PETSC_DEFAULT in 104978bcb3b5SPeter Brune place of an argument. 1050f40b411bSPeter Brune 105178bcb3b5SPeter Brune Level: intermediate 1052516fe3c3SPeter Brune 1053f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances() 1054f40b411bSPeter Brune @*/ 1055f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 10566a388c36SPeter Brune { 10576a388c36SPeter Brune PetscFunctionBegin; 1058f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1059d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1060d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1061d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1062d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,atol,5); 1063d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1064d3952184SSatish Balay PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1065d3952184SSatish Balay 1066d3952184SSatish Balay if (steptol!= PETSC_DEFAULT) { 1067ce94432eSBarry Smith if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 10686a388c36SPeter Brune linesearch->steptol = steptol; 1069d3952184SSatish Balay } 1070d3952184SSatish Balay 1071d3952184SSatish Balay if (maxstep!= PETSC_DEFAULT) { 1072ce94432eSBarry Smith if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1073516fe3c3SPeter Brune linesearch->maxstep = maxstep; 1074d3952184SSatish Balay } 1075d3952184SSatish Balay 1076d3952184SSatish Balay if (rtol != PETSC_DEFAULT) { 1077ce94432eSBarry Smith if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol); 1078516fe3c3SPeter Brune linesearch->rtol = rtol; 1079d3952184SSatish Balay } 1080d3952184SSatish Balay 1081d3952184SSatish Balay if (atol != PETSC_DEFAULT) { 1082ce94432eSBarry Smith if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1083516fe3c3SPeter Brune linesearch->atol = atol; 1084d3952184SSatish Balay } 1085d3952184SSatish Balay 1086d3952184SSatish Balay if (ltol != PETSC_DEFAULT) { 1087ce94432eSBarry Smith if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1088516fe3c3SPeter Brune linesearch->ltol = ltol; 1089d3952184SSatish Balay } 1090d3952184SSatish Balay 1091d3952184SSatish Balay if (max_its != PETSC_DEFAULT) { 1092ce94432eSBarry Smith if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1093516fe3c3SPeter Brune linesearch->max_its = max_its; 1094d3952184SSatish Balay } 10956a388c36SPeter Brune PetscFunctionReturn(0); 10966a388c36SPeter Brune } 10976a388c36SPeter Brune 10986a388c36SPeter Brune #undef __FUNCT__ 1099f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping" 1100f40b411bSPeter Brune /*@ 1101f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter. 1102f40b411bSPeter Brune 1103f40b411bSPeter Brune Input Parameters: 11048e557f58SPeter Brune . linesearch - linesearch context 1105f40b411bSPeter Brune 1106f40b411bSPeter Brune Output Parameters: 11078e557f58SPeter Brune . damping - The damping parameter 1108f40b411bSPeter Brune 110978bcb3b5SPeter Brune Level: advanced 1110f40b411bSPeter Brune 111178bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1112f40b411bSPeter Brune @*/ 1113f40b411bSPeter Brune 1114f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 11156a388c36SPeter Brune { 11166a388c36SPeter Brune PetscFunctionBegin; 1117f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11186a388c36SPeter Brune PetscValidPointer(damping, 2); 11196a388c36SPeter Brune *damping = linesearch->damping; 11206a388c36SPeter Brune PetscFunctionReturn(0); 11216a388c36SPeter Brune } 11226a388c36SPeter Brune 11236a388c36SPeter Brune #undef __FUNCT__ 1124f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping" 1125f40b411bSPeter Brune /*@ 1126f1c6b773SPeter Brune SNESLineSearchSetDamping - Sets the line search damping paramter. 1127f40b411bSPeter Brune 1128f40b411bSPeter Brune Input Parameters: 112978bcb3b5SPeter Brune . linesearch - linesearch context 113078bcb3b5SPeter Brune . damping - The damping parameter 1131f40b411bSPeter Brune 1132f40b411bSPeter Brune Level: intermediate 1133f40b411bSPeter Brune 1134cd7522eaSPeter Brune Notes: 1135cd7522eaSPeter Brune The basic line search merely takes the update step scaled by the damping parameter. 1136cd7522eaSPeter Brune The use of the damping parameter in the l2 and cp line searches is much more subtle; 113778bcb3b5SPeter Brune it is used as a starting point in calculating the secant step. However, the eventual 1138cd7522eaSPeter Brune step may be of greater length than the damping parameter. In the bt line search it is 1139cd7522eaSPeter Brune used as the maximum possible step length, as the bt line search only backtracks. 1140cd7522eaSPeter Brune 1141f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping() 1142f40b411bSPeter Brune @*/ 1143f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 11446a388c36SPeter Brune { 11456a388c36SPeter Brune PetscFunctionBegin; 1146f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 11476a388c36SPeter Brune linesearch->damping = damping; 11486a388c36SPeter Brune PetscFunctionReturn(0); 11496a388c36SPeter Brune } 11506a388c36SPeter Brune 11516a388c36SPeter Brune #undef __FUNCT__ 115259405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder" 115359405d5eSPeter Brune /*@ 115459405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order. 115559405d5eSPeter Brune 115659405d5eSPeter Brune Input Parameters: 115778bcb3b5SPeter Brune . linesearch - linesearch context 115859405d5eSPeter Brune 115959405d5eSPeter Brune Output Parameters: 11608e557f58SPeter Brune . order - The order 116159405d5eSPeter Brune 116278bcb3b5SPeter Brune Possible Values for order: 11633c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 11643c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 11653c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 116678bcb3b5SPeter Brune 116759405d5eSPeter Brune Level: intermediate 116859405d5eSPeter Brune 116959405d5eSPeter Brune .seealso: SNESLineSearchSetOrder() 117059405d5eSPeter Brune @*/ 117159405d5eSPeter Brune 1172b000cd8dSPeter Brune PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 117359405d5eSPeter Brune { 117459405d5eSPeter Brune PetscFunctionBegin; 117559405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 117659405d5eSPeter Brune PetscValidPointer(order, 2); 117759405d5eSPeter Brune *order = linesearch->order; 117859405d5eSPeter Brune PetscFunctionReturn(0); 117959405d5eSPeter Brune } 118059405d5eSPeter Brune 118159405d5eSPeter Brune #undef __FUNCT__ 118259405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder" 118359405d5eSPeter Brune /*@ 118459405d5eSPeter Brune SNESLineSearchSetOrder - Sets the line search damping paramter. 118559405d5eSPeter Brune 118659405d5eSPeter Brune Input Parameters: 118778bcb3b5SPeter Brune . linesearch - linesearch context 118878bcb3b5SPeter Brune . order - The damping parameter 118959405d5eSPeter Brune 119059405d5eSPeter Brune Level: intermediate 119159405d5eSPeter Brune 119278bcb3b5SPeter Brune Possible Values for order: 11933c7d6663SPeter Brune + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 11943c7d6663SPeter Brune . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 11953c7d6663SPeter Brune - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 119678bcb3b5SPeter Brune 119759405d5eSPeter Brune Notes: 119859405d5eSPeter Brune Variable orders are supported by the following line searches: 119978bcb3b5SPeter Brune + bt - cubic and quadratic 120078bcb3b5SPeter Brune - cp - linear and quadratic 120159405d5eSPeter Brune 120259405d5eSPeter Brune .seealso: SNESLineSearchGetOrder() 120359405d5eSPeter Brune @*/ 1204b000cd8dSPeter Brune PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 120559405d5eSPeter Brune { 120659405d5eSPeter Brune PetscFunctionBegin; 120759405d5eSPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 120859405d5eSPeter Brune linesearch->order = order; 120959405d5eSPeter Brune PetscFunctionReturn(0); 121059405d5eSPeter Brune } 121159405d5eSPeter Brune 121259405d5eSPeter Brune #undef __FUNCT__ 1213f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms" 1214f40b411bSPeter Brune /*@ 1215f1c6b773SPeter Brune SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1216f40b411bSPeter Brune 1217f40b411bSPeter Brune Input Parameters: 121878bcb3b5SPeter Brune . linesearch - linesearch context 1219f40b411bSPeter Brune 1220f40b411bSPeter Brune Output Parameters: 1221f40b411bSPeter Brune + xnorm - The norm of the current solution 1222f40b411bSPeter Brune . fnorm - The norm of the current function 1223f40b411bSPeter Brune - ynorm - The norm of the current update 1224f40b411bSPeter Brune 1225cd7522eaSPeter Brune Notes: 1226cd7522eaSPeter Brune This function is mainly called from SNES implementations. 1227cd7522eaSPeter Brune 122878bcb3b5SPeter Brune Level: developer 1229f40b411bSPeter Brune 1230f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1231f40b411bSPeter Brune @*/ 1232f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1233bf7f4e0aSPeter Brune { 1234bf7f4e0aSPeter Brune PetscFunctionBegin; 1235f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1236f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm; 1237f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm; 1238f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm; 1239bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1240bf7f4e0aSPeter Brune } 1241bf7f4e0aSPeter Brune 1242e7058c64SPeter Brune #undef __FUNCT__ 1243f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms" 1244f40b411bSPeter Brune /*@ 1245f1c6b773SPeter Brune SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1246f40b411bSPeter Brune 1247f40b411bSPeter Brune Input Parameters: 124878bcb3b5SPeter Brune + linesearch - linesearch context 1249f40b411bSPeter Brune . xnorm - The norm of the current solution 1250f40b411bSPeter Brune . fnorm - The norm of the current function 1251f40b411bSPeter Brune - ynorm - The norm of the current update 1252f40b411bSPeter Brune 125378bcb3b5SPeter Brune Level: advanced 1254f40b411bSPeter Brune 1255f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1256f40b411bSPeter Brune @*/ 1257f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 12586a388c36SPeter Brune { 12596a388c36SPeter Brune PetscFunctionBegin; 1260f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 12616a388c36SPeter Brune linesearch->xnorm = xnorm; 12626a388c36SPeter Brune linesearch->fnorm = fnorm; 12636a388c36SPeter Brune linesearch->ynorm = ynorm; 12646a388c36SPeter Brune PetscFunctionReturn(0); 12656a388c36SPeter Brune } 12666a388c36SPeter Brune 12676a388c36SPeter Brune #undef __FUNCT__ 1268f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms" 1269f40b411bSPeter Brune /*@ 1270f1c6b773SPeter Brune SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1271f40b411bSPeter Brune 1272f40b411bSPeter Brune Input Parameters: 127378bcb3b5SPeter Brune . linesearch - linesearch context 1274f40b411bSPeter Brune 1275f40b411bSPeter Brune Options Database Keys: 12768e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 1277f40b411bSPeter Brune 1278f40b411bSPeter Brune Level: intermediate 1279f40b411bSPeter Brune 128078bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1281f40b411bSPeter Brune @*/ 1282f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 12836a388c36SPeter Brune { 12846a388c36SPeter Brune PetscErrorCode ierr; 12859bd66eb0SPeter Brune SNES snes; 1286bf388a1fSBarry Smith 12876a388c36SPeter Brune PetscFunctionBegin; 12886a388c36SPeter Brune if (linesearch->norms) { 12899bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1290f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 12919bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12929bd66eb0SPeter Brune ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12939bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 12949bd66eb0SPeter Brune } else { 12956a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12966a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 12976a388c36SPeter Brune ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 12986a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 12996a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 13006a388c36SPeter Brune ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 13016a388c36SPeter Brune } 13029bd66eb0SPeter Brune } 13036a388c36SPeter Brune PetscFunctionReturn(0); 13046a388c36SPeter Brune } 13056a388c36SPeter Brune 13066f263ca3SPeter Brune #undef __FUNCT__ 13076f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms" 13086f263ca3SPeter Brune /*@ 13096f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 13106f263ca3SPeter Brune 13116f263ca3SPeter Brune Input Parameters: 131278bcb3b5SPeter Brune + linesearch - linesearch context 131378bcb3b5SPeter Brune - flg - indicates whether or not to compute norms 13146f263ca3SPeter Brune 13156f263ca3SPeter Brune Options Database Keys: 13168e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off 13176f263ca3SPeter Brune 13186f263ca3SPeter Brune Notes: 13191a4f838cSPeter Brune This is most relevant to the SNESLINESEARCHBASIC line search type. 13206f263ca3SPeter Brune 13216f263ca3SPeter Brune Level: intermediate 13226f263ca3SPeter Brune 13231a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 13246f263ca3SPeter Brune @*/ 13256f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 13266f263ca3SPeter Brune { 13276f263ca3SPeter Brune PetscFunctionBegin; 13286f263ca3SPeter Brune linesearch->norms = flg; 13296f263ca3SPeter Brune PetscFunctionReturn(0); 13306f263ca3SPeter Brune } 13316f263ca3SPeter Brune 13326a388c36SPeter Brune #undef __FUNCT__ 1333f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs" 1334f40b411bSPeter Brune /*@ 1335f1c6b773SPeter Brune SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1336f40b411bSPeter Brune 1337f40b411bSPeter Brune Input Parameters: 133878bcb3b5SPeter Brune . linesearch - linesearch context 1339f40b411bSPeter Brune 1340f40b411bSPeter Brune Output Parameters: 13416232e825SPeter Brune + X - Solution vector 13426232e825SPeter Brune . F - Function vector 13436232e825SPeter Brune . Y - Search direction vector 13446232e825SPeter Brune . W - Solution work vector 13456232e825SPeter Brune - G - Function work vector 13466232e825SPeter Brune 13477bba9028SPeter Brune Notes: 13487bba9028SPeter Brune At the beginning of a line search application, X should contain a 13496232e825SPeter Brune solution and the vector F the function computed at X. At the end of the 13506232e825SPeter Brune line search application, X should contain the new solution, and F the 13516232e825SPeter Brune function evaluated at the new solution. 1352f40b411bSPeter Brune 13532a7a6963SBarry Smith These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 13542a7a6963SBarry Smith 135578bcb3b5SPeter Brune Level: advanced 1356f40b411bSPeter Brune 1357f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1358f40b411bSPeter Brune @*/ 1359bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1360bf388a1fSBarry Smith { 13616a388c36SPeter Brune PetscFunctionBegin; 1362f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 13636a388c36SPeter Brune if (X) { 13646a388c36SPeter Brune PetscValidPointer(X, 2); 13656a388c36SPeter Brune *X = linesearch->vec_sol; 13666a388c36SPeter Brune } 13676a388c36SPeter Brune if (F) { 13686a388c36SPeter Brune PetscValidPointer(F, 3); 13696a388c36SPeter Brune *F = linesearch->vec_func; 13706a388c36SPeter Brune } 13716a388c36SPeter Brune if (Y) { 13726a388c36SPeter Brune PetscValidPointer(Y, 4); 13736a388c36SPeter Brune *Y = linesearch->vec_update; 13746a388c36SPeter Brune } 13756a388c36SPeter Brune if (W) { 13766a388c36SPeter Brune PetscValidPointer(W, 5); 13776a388c36SPeter Brune *W = linesearch->vec_sol_new; 13786a388c36SPeter Brune } 13796a388c36SPeter Brune if (G) { 13806a388c36SPeter Brune PetscValidPointer(G, 6); 13816a388c36SPeter Brune *G = linesearch->vec_func_new; 13826a388c36SPeter Brune } 13836a388c36SPeter Brune PetscFunctionReturn(0); 13846a388c36SPeter Brune } 13856a388c36SPeter Brune 13866a388c36SPeter Brune #undef __FUNCT__ 1387f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs" 1388f40b411bSPeter Brune /*@ 1389f1c6b773SPeter Brune SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1390f40b411bSPeter Brune 1391f40b411bSPeter Brune Input Parameters: 139278bcb3b5SPeter Brune + linesearch - linesearch context 13936232e825SPeter Brune . X - Solution vector 13946232e825SPeter Brune . F - Function vector 13956232e825SPeter Brune . Y - Search direction vector 13966232e825SPeter Brune . W - Solution work vector 13976232e825SPeter Brune - G - Function work vector 1398f40b411bSPeter Brune 139978bcb3b5SPeter Brune Level: advanced 1400f40b411bSPeter Brune 1401f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1402f40b411bSPeter Brune @*/ 1403bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1404bf388a1fSBarry Smith { 14056a388c36SPeter Brune PetscFunctionBegin; 1406f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 14076a388c36SPeter Brune if (X) { 14086a388c36SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 14096a388c36SPeter Brune linesearch->vec_sol = X; 14106a388c36SPeter Brune } 14116a388c36SPeter Brune if (F) { 14126a388c36SPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 14136a388c36SPeter Brune linesearch->vec_func = F; 14146a388c36SPeter Brune } 14156a388c36SPeter Brune if (Y) { 14166a388c36SPeter Brune PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 14176a388c36SPeter Brune linesearch->vec_update = Y; 14186a388c36SPeter Brune } 14196a388c36SPeter Brune if (W) { 14206a388c36SPeter Brune PetscValidHeaderSpecific(W,VEC_CLASSID,5); 14216a388c36SPeter Brune linesearch->vec_sol_new = W; 14226a388c36SPeter Brune } 14236a388c36SPeter Brune if (G) { 14246a388c36SPeter Brune PetscValidHeaderSpecific(G,VEC_CLASSID,6); 14256a388c36SPeter Brune linesearch->vec_func_new = G; 14266a388c36SPeter Brune } 14276a388c36SPeter Brune PetscFunctionReturn(0); 14286a388c36SPeter Brune } 14296a388c36SPeter Brune 14306a388c36SPeter Brune #undef __FUNCT__ 1431f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1432e7058c64SPeter Brune /*@C 1433f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1434e7058c64SPeter Brune SNES options in the database. 1435e7058c64SPeter Brune 1436cd7522eaSPeter Brune Logically Collective on SNESLineSearch 1437e7058c64SPeter Brune 1438e7058c64SPeter Brune Input Parameters: 1439e7058c64SPeter Brune + snes - the SNES context 1440e7058c64SPeter Brune - prefix - the prefix to prepend to all option names 1441e7058c64SPeter Brune 1442e7058c64SPeter Brune Notes: 1443e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name. 1444e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen. 1445e7058c64SPeter Brune 1446e7058c64SPeter Brune Level: advanced 1447e7058c64SPeter Brune 1448f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database 1449e7058c64SPeter Brune 1450e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix() 1451e7058c64SPeter Brune @*/ 1452f1c6b773SPeter Brune PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1453e7058c64SPeter Brune { 1454e7058c64SPeter Brune PetscErrorCode ierr; 1455e7058c64SPeter Brune 1456e7058c64SPeter Brune PetscFunctionBegin; 1457f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1458e7058c64SPeter Brune ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1459e7058c64SPeter Brune PetscFunctionReturn(0); 1460e7058c64SPeter Brune } 1461e7058c64SPeter Brune 1462e7058c64SPeter Brune #undef __FUNCT__ 1463f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1464e7058c64SPeter Brune /*@C 1465f1c6b773SPeter Brune SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1466f1c6b773SPeter Brune SNESLineSearch options in the database. 1467e7058c64SPeter Brune 1468e7058c64SPeter Brune Not Collective 1469e7058c64SPeter Brune 1470e7058c64SPeter Brune Input Parameter: 1471cd7522eaSPeter Brune . linesearch - the SNESLineSearch context 1472e7058c64SPeter Brune 1473e7058c64SPeter Brune Output Parameter: 1474e7058c64SPeter Brune . prefix - pointer to the prefix string used 1475e7058c64SPeter Brune 14768e557f58SPeter Brune Notes: 14778e557f58SPeter Brune On the fortran side, the user should pass in a string 'prefix' of 1478e7058c64SPeter Brune sufficient length to hold the prefix. 1479e7058c64SPeter Brune 1480e7058c64SPeter Brune Level: advanced 1481e7058c64SPeter Brune 1482f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database 1483e7058c64SPeter Brune 1484e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix() 1485e7058c64SPeter Brune @*/ 1486f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1487e7058c64SPeter Brune { 1488e7058c64SPeter Brune PetscErrorCode ierr; 1489e7058c64SPeter Brune 1490e7058c64SPeter Brune PetscFunctionBegin; 1491f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1492e7058c64SPeter Brune ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1493e7058c64SPeter Brune PetscFunctionReturn(0); 1494e7058c64SPeter Brune } 1495bf7f4e0aSPeter Brune 1496bf7f4e0aSPeter Brune #undef __FUNCT__ 1497fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs" 14988d359177SBarry Smith /*@C 1499fa0ddf94SBarry Smith SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1500f40b411bSPeter Brune 1501f40b411bSPeter Brune Input Parameter: 1502f1c6b773SPeter Brune + linesearch - the SNESLineSearch context 1503f40b411bSPeter Brune - nwork - the number of work vectors 1504f40b411bSPeter Brune 1505f40b411bSPeter Brune Level: developer 1506f40b411bSPeter Brune 1507fa0ddf94SBarry Smith Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1508cd7522eaSPeter Brune 1509f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector 1510f40b411bSPeter Brune 1511fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs() 1512f40b411bSPeter Brune @*/ 1513fa0ddf94SBarry Smith PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1514bf7f4e0aSPeter Brune { 1515bf7f4e0aSPeter Brune PetscErrorCode ierr; 1516bf388a1fSBarry Smith 1517bf7f4e0aSPeter Brune PetscFunctionBegin; 1518bf7f4e0aSPeter Brune if (linesearch->vec_sol) { 1519bf7f4e0aSPeter Brune ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 15208d359177SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1521bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1522bf7f4e0aSPeter Brune } 1523bf7f4e0aSPeter Brune 1524bf7f4e0aSPeter Brune #undef __FUNCT__ 1525422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchGetReason" 1526f40b411bSPeter Brune /*@ 1527422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1528f40b411bSPeter Brune 1529f40b411bSPeter Brune Input Parameters: 153078bcb3b5SPeter Brune . linesearch - linesearch context 1531f40b411bSPeter Brune 1532f40b411bSPeter Brune Output Parameters: 1533422a814eSBarry Smith . result - The success or failure status 1534f40b411bSPeter Brune 1535cd7522eaSPeter Brune Notes: 1536c98378a5SBarry Smith This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1537cd7522eaSPeter Brune (and set the SNES convergence accordingly). 1538cd7522eaSPeter Brune 1539f40b411bSPeter Brune Level: intermediate 1540f40b411bSPeter Brune 1541422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason 1542f40b411bSPeter Brune @*/ 1543422a814eSBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1544bf7f4e0aSPeter Brune { 1545bf7f4e0aSPeter Brune PetscFunctionBegin; 1546f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1547422a814eSBarry Smith PetscValidPointer(result, 2); 1548422a814eSBarry Smith *result = linesearch->result; 1549bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1550bf7f4e0aSPeter Brune } 1551bf7f4e0aSPeter Brune 1552bf7f4e0aSPeter Brune #undef __FUNCT__ 1553422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchSetReason" 1554f40b411bSPeter Brune /*@ 1555422a814eSBarry Smith SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1556f40b411bSPeter Brune 1557f40b411bSPeter Brune Input Parameters: 155878bcb3b5SPeter Brune + linesearch - linesearch context 1559422a814eSBarry Smith - result - The success or failure status 1560f40b411bSPeter Brune 1561cd7522eaSPeter Brune Notes: 1562cd7522eaSPeter Brune This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1563cd7522eaSPeter Brune the success or failure of the line search method. 1564cd7522eaSPeter Brune 156578bcb3b5SPeter Brune Level: developer 1566f40b411bSPeter Brune 1567422a814eSBarry Smith .seealso: SNESLineSearchGetSResult() 1568f40b411bSPeter Brune @*/ 1569422a814eSBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 15706a388c36SPeter Brune { 15716a388c36SPeter Brune PetscFunctionBegin; 15725fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1573422a814eSBarry Smith linesearch->result = result; 15746a388c36SPeter Brune PetscFunctionReturn(0); 15756a388c36SPeter Brune } 15766a388c36SPeter Brune 15776a388c36SPeter Brune #undef __FUNCT__ 1578f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions" 15799bd66eb0SPeter Brune /*@C 1580f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 15819bd66eb0SPeter Brune 15829bd66eb0SPeter Brune Input Parameters: 15839bd66eb0SPeter Brune + snes - nonlinear context obtained from SNESCreate() 15849bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds 15859bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 15869bd66eb0SPeter Brune 15879bd66eb0SPeter Brune Logically Collective on SNES 15889bd66eb0SPeter Brune 15899bd66eb0SPeter Brune Calling sequence of projectfunc: 15909bd66eb0SPeter Brune .vb 15919bd66eb0SPeter Brune projectfunc (SNES snes, Vec X) 15929bd66eb0SPeter Brune .ve 15939bd66eb0SPeter Brune 15949bd66eb0SPeter Brune Input parameters for projectfunc: 15959bd66eb0SPeter Brune + snes - nonlinear context 15969bd66eb0SPeter Brune - X - current solution 15979bd66eb0SPeter Brune 1598cd7522eaSPeter Brune Output parameters for projectfunc: 15999bd66eb0SPeter Brune . X - Projected solution 16009bd66eb0SPeter Brune 16019bd66eb0SPeter Brune Calling sequence of normfunc: 16029bd66eb0SPeter Brune .vb 16039bd66eb0SPeter Brune projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 16049bd66eb0SPeter Brune .ve 16059bd66eb0SPeter Brune 1606cd7522eaSPeter Brune Input parameters for normfunc: 16079bd66eb0SPeter Brune + snes - nonlinear context 16089bd66eb0SPeter Brune . X - current solution 16099bd66eb0SPeter Brune - F - current residual 16109bd66eb0SPeter Brune 1611cd7522eaSPeter Brune Output parameters for normfunc: 16129bd66eb0SPeter Brune . fnorm - VI-specific norm of the function 16139bd66eb0SPeter Brune 1614cd7522eaSPeter Brune Notes: 1615cd7522eaSPeter Brune The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1616cd7522eaSPeter Brune 1617cd7522eaSPeter Brune The VI solvers require special evaluation of the function norm such that the norm is only calculated 1618cd7522eaSPeter Brune on the inactive set. This should be implemented by normfunc. 16199bd66eb0SPeter Brune 16209bd66eb0SPeter Brune Level: developer 16219bd66eb0SPeter Brune 16229bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search 16239bd66eb0SPeter Brune 1624f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 16259bd66eb0SPeter Brune @*/ 1626f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 16279bd66eb0SPeter Brune { 16289bd66eb0SPeter Brune PetscFunctionBegin; 1629f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 16309bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc; 16319bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc; 16329bd66eb0SPeter Brune PetscFunctionReturn(0); 16339bd66eb0SPeter Brune } 16349bd66eb0SPeter Brune 16359bd66eb0SPeter Brune #undef __FUNCT__ 1636f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions" 16379bd66eb0SPeter Brune /*@C 1638f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 16399bd66eb0SPeter Brune 16409bd66eb0SPeter Brune Input Parameters: 1641907376e6SBarry Smith . linesearch - the line search context, obtain with SNESGetLineSearch() 16429bd66eb0SPeter Brune 16439bd66eb0SPeter Brune Output Parameters: 16449bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds 16459bd66eb0SPeter Brune - normfunc - function for computing the norm of an active set 16469bd66eb0SPeter Brune 16479bd66eb0SPeter Brune Logically Collective on SNES 16489bd66eb0SPeter Brune 16499bd66eb0SPeter Brune Level: developer 16509bd66eb0SPeter Brune 16519bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search 16529bd66eb0SPeter Brune 1653f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 16549bd66eb0SPeter Brune @*/ 1655f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 16569bd66eb0SPeter Brune { 16579bd66eb0SPeter Brune PetscFunctionBegin; 16589bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject; 16599bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm; 16609bd66eb0SPeter Brune PetscFunctionReturn(0); 16619bd66eb0SPeter Brune } 16629bd66eb0SPeter Brune 16639bd66eb0SPeter Brune #undef __FUNCT__ 1664f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister" 1665bf7f4e0aSPeter Brune /*@C 16661c84c290SBarry Smith SNESLineSearchRegister - See SNESLineSearchRegister() 1667bf7f4e0aSPeter Brune 1668bf7f4e0aSPeter Brune Level: advanced 1669bf7f4e0aSPeter Brune @*/ 1670bdf89e91SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1671bf7f4e0aSPeter Brune { 1672bf7f4e0aSPeter Brune PetscErrorCode ierr; 1673bf7f4e0aSPeter Brune 1674bf7f4e0aSPeter Brune PetscFunctionBegin; 1675a240a19fSJed Brown ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr); 1676bf7f4e0aSPeter Brune PetscFunctionReturn(0); 1677bf7f4e0aSPeter Brune } 1678