xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 86d74e61e44fcf73060f77fc8677142765f9cd55)
19e764e56SPeter Brune #include <private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
36188f407SPeter Brune PetscBool  PetscLineSearchRegisterAllCalled = PETSC_FALSE;
46188f407SPeter Brune PetscFList PetscLineSearchList              = PETSC_NULL;
5bf7f4e0aSPeter Brune 
66188f407SPeter Brune PetscClassId   PETSCLINESEARCH_CLASSID;
76188f407SPeter Brune PetscLogEvent  PetscLineSearch_Apply;
8bf7f4e0aSPeter Brune 
9bf7f4e0aSPeter Brune #undef __FUNCT__
106188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate"
11f40b411bSPeter Brune /*@
126188f407SPeter Brune    PetscLineSearchCreate - Creates the line search.
13f40b411bSPeter Brune 
14f40b411bSPeter Brune    Collective on LineSearch
15f40b411bSPeter Brune 
16f40b411bSPeter Brune    Input Parameters:
17f40b411bSPeter Brune .  comm - MPI communicator for the line search
18f40b411bSPeter Brune 
19f40b411bSPeter Brune    Output Parameters:
20f40b411bSPeter Brune .  outlinesearch - the line search instance.
21f40b411bSPeter Brune 
22f40b411bSPeter Brune    Level: Beginner
23f40b411bSPeter Brune 
24f40b411bSPeter Brune    .keywords: LineSearch, Create
25f40b411bSPeter Brune 
26f40b411bSPeter Brune    .seealso: LineSearchDestroy()
27f40b411bSPeter Brune @*/
28f40b411bSPeter Brune 
296188f407SPeter Brune PetscErrorCode PetscLineSearchCreate(MPI_Comm comm, PetscLineSearch * outlinesearch) {
30bf7f4e0aSPeter Brune   PetscErrorCode ierr;
316188f407SPeter Brune   PetscLineSearch     linesearch;
32bf7f4e0aSPeter Brune   PetscFunctionBegin;
336188f407SPeter Brune   ierr = PetscHeaderCreate(linesearch, _p_LineSearch,struct _LineSearchOps,PETSCLINESEARCH_CLASSID, 0,
346188f407SPeter Brune                            "LineSearch","Line-search method","LineSearch",comm,PetscLineSearchDestroy,PetscLineSearchView);CHKERRQ(ierr);
35bf7f4e0aSPeter Brune 
36bf7f4e0aSPeter Brune   linesearch->ops->precheckstep = PETSC_NULL;
37bf7f4e0aSPeter Brune   linesearch->ops->postcheckstep = PETSC_NULL;
38bf7f4e0aSPeter Brune 
39bf7f4e0aSPeter Brune   linesearch->lambda        = 1.0;
40bf7f4e0aSPeter Brune   linesearch->fnorm         = 1.0;
41bf7f4e0aSPeter Brune   linesearch->ynorm         = 1.0;
42bf7f4e0aSPeter Brune   linesearch->xnorm         = 1.0;
43bf7f4e0aSPeter Brune   linesearch->success       = PETSC_TRUE;
44bf7f4e0aSPeter Brune   linesearch->norms         = PETSC_TRUE;
45bf7f4e0aSPeter Brune   linesearch->keeplambda    = PETSC_FALSE;
46bf7f4e0aSPeter Brune   linesearch->damping       = 1.0;
47bf7f4e0aSPeter Brune   linesearch->maxstep       = 1e8;
48bf7f4e0aSPeter Brune   linesearch->steptol       = 1e-12;
49516fe3c3SPeter Brune   linesearch->rtol          = 1e-8;
50516fe3c3SPeter Brune   linesearch->atol          = 1e-15;
51516fe3c3SPeter Brune   linesearch->ltol          = 1e-8;
52bf7f4e0aSPeter Brune   linesearch->precheckctx   = PETSC_NULL;
53bf7f4e0aSPeter Brune   linesearch->postcheckctx  = PETSC_NULL;
54516fe3c3SPeter Brune   linesearch->max_its       = 3;
55bf7f4e0aSPeter Brune   linesearch->setupcalled   = PETSC_FALSE;
56bf7f4e0aSPeter Brune   *outlinesearch            = linesearch;
57bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
58bf7f4e0aSPeter Brune }
59bf7f4e0aSPeter Brune 
60bf7f4e0aSPeter Brune #undef __FUNCT__
616188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetUp"
62f40b411bSPeter Brune /*@
636188f407SPeter Brune    PetscLineSearchSetUp - Prepares the line search for being applied.
64f40b411bSPeter Brune 
65f40b411bSPeter Brune    Collective on LineSearch
66f40b411bSPeter Brune 
67f40b411bSPeter Brune    Input Parameters:
68f40b411bSPeter Brune .  linesearch - The LineSearch instance.
69f40b411bSPeter Brune 
70f40b411bSPeter Brune    Level: Intermediate
71f40b411bSPeter Brune 
726188f407SPeter Brune    .keywords: PetscLineSearch, SetUp
73f40b411bSPeter Brune 
746188f407SPeter Brune    .seealso: PetscLineSearchReset()
75f40b411bSPeter Brune @*/
76f40b411bSPeter Brune 
776188f407SPeter Brune PetscErrorCode PetscLineSearchSetUp(PetscLineSearch linesearch) {
78bf7f4e0aSPeter Brune   PetscErrorCode ierr;
79bf7f4e0aSPeter Brune   PetscFunctionBegin;
80bf7f4e0aSPeter Brune 
81bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
826188f407SPeter Brune     ierr = PetscLineSearchSetType(linesearch,PETSCLINESEARCHBASIC);CHKERRQ(ierr);
83bf7f4e0aSPeter Brune   }
84bf7f4e0aSPeter Brune 
85bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
86bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
87bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_func, &linesearch->vec_func_new);CHKERRQ(ierr);
88bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
89bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
90bf7f4e0aSPeter Brune     }
91bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
92bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
93bf7f4e0aSPeter Brune   }
94bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
95bf7f4e0aSPeter Brune }
96bf7f4e0aSPeter Brune 
97bf7f4e0aSPeter Brune #undef __FUNCT__
986188f407SPeter Brune #define __FUNCT__ "PetscLineSearchReset"
99f40b411bSPeter Brune 
100f40b411bSPeter Brune /*@
1016188f407SPeter Brune    PetscLineSearchReset - Tears down the structures required for application
102f40b411bSPeter Brune 
1036188f407SPeter Brune    Collective on PetscLineSearch
104f40b411bSPeter Brune 
105f40b411bSPeter Brune    Input Parameters:
106f40b411bSPeter Brune .  linesearch - The LineSearch instance.
107f40b411bSPeter Brune 
108f40b411bSPeter Brune    Level: Intermediate
109f40b411bSPeter Brune 
1106188f407SPeter Brune    .keywords: PetscLineSearch, Create
111f40b411bSPeter Brune 
1126188f407SPeter Brune    .seealso: PetscLineSearchSetUp()
113f40b411bSPeter Brune @*/
114f40b411bSPeter Brune 
1156188f407SPeter Brune PetscErrorCode PetscLineSearchReset(PetscLineSearch linesearch) {
116bf7f4e0aSPeter Brune   PetscErrorCode ierr;
117bf7f4e0aSPeter Brune   PetscFunctionBegin;
118bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
119bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
120bf7f4e0aSPeter Brune   }
121bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
122bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
123bf7f4e0aSPeter Brune 
124bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
125bf7f4e0aSPeter Brune   linesearch->nwork = 0;
126bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
127bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
128bf7f4e0aSPeter Brune }
129bf7f4e0aSPeter Brune 
130*86d74e61SPeter Brune 
131*86d74e61SPeter Brune #undef __FUNCT__
132*86d74e61SPeter Brune #define __FUNCT__ "PetscLineSearchSetPreCheck"
133*86d74e61SPeter Brune /*@C
134*86d74e61SPeter Brune    PetscLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
135*86d74e61SPeter Brune 
136*86d74e61SPeter Brune    Logically Collective on PetscLineSearch
137*86d74e61SPeter Brune 
138*86d74e61SPeter Brune    Input Parameters:
139*86d74e61SPeter Brune +  linesearch - the PetscLineSearch context
140*86d74e61SPeter Brune .  func       - [optional] function evaluation routine
141*86d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
142*86d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
143*86d74e61SPeter Brune 
144*86d74e61SPeter Brune    Calling sequence of func:
145*86d74e61SPeter Brune $    func (PetscLineSearch snes,Vec x,Vec y, PetscBool *changed);
146*86d74e61SPeter Brune 
147*86d74e61SPeter Brune +  x - solution vector
148*86d74e61SPeter Brune .  y - search direction vector
149*86d74e61SPeter Brune -  changed - flag to indicate the precheck changed something.
150*86d74e61SPeter Brune 
151*86d74e61SPeter Brune    Level: intermediate
152*86d74e61SPeter Brune 
153*86d74e61SPeter Brune .keywords: set, linesearch, pre-check
154*86d74e61SPeter Brune 
155*86d74e61SPeter Brune .seealso: PetscLineSearchSetPostCheck()
156*86d74e61SPeter Brune @*/
157*86d74e61SPeter Brune PetscErrorCode  PetscLineSearchSetPreCheck(PetscLineSearch linesearch, PetscLineSearchPreCheckFunc func,void *ctx)
158*86d74e61SPeter Brune {
159*86d74e61SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
160*86d74e61SPeter Brune   if (func) linesearch->ops->precheckstep = func;
161*86d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
162*86d74e61SPeter Brune   PetscFunctionReturn(0);
163*86d74e61SPeter Brune }
164*86d74e61SPeter Brune 
165*86d74e61SPeter Brune 
166*86d74e61SPeter Brune #undef __FUNCT__
167*86d74e61SPeter Brune #define __FUNCT__ "PetscLineSearchGetPreCheck"
168*86d74e61SPeter Brune /*@C
169*86d74e61SPeter Brune    PetscLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
170*86d74e61SPeter Brune 
171*86d74e61SPeter Brune    Input Parameters:
172*86d74e61SPeter Brune .  linesearch - the PetscLineSearch context
173*86d74e61SPeter Brune 
174*86d74e61SPeter Brune    Output Parameters:
175*86d74e61SPeter Brune +  func       - [optional] function evaluation routine
176*86d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
177*86d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
178*86d74e61SPeter Brune 
179*86d74e61SPeter Brune    Level: intermediate
180*86d74e61SPeter Brune 
181*86d74e61SPeter Brune .keywords: get, linesearch, pre-check
182*86d74e61SPeter Brune 
183*86d74e61SPeter Brune .seealso: PetscLineSearchGetPostCheck(), PetscLineSearchSetPreCheck()
184*86d74e61SPeter Brune @*/
185*86d74e61SPeter Brune PetscErrorCode  PetscLineSearchGetPreCheck(PetscLineSearch linesearch, PetscLineSearchPreCheckFunc *func,void **ctx)
186*86d74e61SPeter Brune {
187*86d74e61SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
188*86d74e61SPeter Brune   if (func) *func = linesearch->ops->precheckstep;
189*86d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
190*86d74e61SPeter Brune   PetscFunctionReturn(0);
191*86d74e61SPeter Brune }
192*86d74e61SPeter Brune 
193*86d74e61SPeter Brune 
194*86d74e61SPeter Brune #undef __FUNCT__
195*86d74e61SPeter Brune #define __FUNCT__ "PetscLineSearchSetPostCheck"
196*86d74e61SPeter Brune /*@C
197*86d74e61SPeter Brune    PetscLineSearchSetPostCheck - Sets a post-check function for the line search routine.
198*86d74e61SPeter Brune 
199*86d74e61SPeter Brune    Logically Collective on PetscLineSearch
200*86d74e61SPeter Brune 
201*86d74e61SPeter Brune    Input Parameters:
202*86d74e61SPeter Brune +  linesearch - the PetscLineSearch context
203*86d74e61SPeter Brune .  func       - [optional] function evaluation routine
204*86d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
205*86d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
206*86d74e61SPeter Brune 
207*86d74e61SPeter Brune    Calling sequence of func:
208*86d74e61SPeter Brune $    func (PetscLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y);
209*86d74e61SPeter Brune 
210*86d74e61SPeter Brune +  x - old solution vector
211*86d74e61SPeter Brune .  y - search direction vector
212*86d74e61SPeter Brune .  w - new solution vector
213*86d74e61SPeter Brune .  changed_y - indicates that the line search changed y.
214*86d74e61SPeter Brune .  changed_w - indicates that the line search changed w.
215*86d74e61SPeter Brune 
216*86d74e61SPeter Brune    Level: intermediate
217*86d74e61SPeter Brune 
218*86d74e61SPeter Brune .keywords: set, linesearch, post-check
219*86d74e61SPeter Brune 
220*86d74e61SPeter Brune .seealso: PetscLineSearchSetPreCheck()
221*86d74e61SPeter Brune @*/
222*86d74e61SPeter Brune PetscErrorCode  PetscLineSearchSetPostCheck(PetscLineSearch linesearch, PetscLineSearchPostCheckFunc func,void *ctx)
223*86d74e61SPeter Brune {
224*86d74e61SPeter Brune   PetscFunctionBegin;
225*86d74e61SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
226*86d74e61SPeter Brune   if (func) linesearch->ops->postcheckstep = func;
227*86d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
228*86d74e61SPeter Brune   PetscFunctionReturn(0);
229*86d74e61SPeter Brune }
230*86d74e61SPeter Brune 
231*86d74e61SPeter Brune 
232*86d74e61SPeter Brune #undef __FUNCT__
233*86d74e61SPeter Brune #define __FUNCT__ "PetscLineSearchGetPostCheck"
234*86d74e61SPeter Brune /*@C
235*86d74e61SPeter Brune    PetscLineSearchGetPostCheck - Gets the post-check function for the line search routine.
236*86d74e61SPeter Brune 
237*86d74e61SPeter Brune    Input Parameters:
238*86d74e61SPeter Brune .  linesearch - the PetscLineSearch context
239*86d74e61SPeter Brune 
240*86d74e61SPeter Brune    Output Parameters:
241*86d74e61SPeter Brune +  func       - [optional] function evaluation routine
242*86d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
243*86d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
244*86d74e61SPeter Brune 
245*86d74e61SPeter Brune    Level: intermediate
246*86d74e61SPeter Brune 
247*86d74e61SPeter Brune .keywords: get, linesearch, post-check
248*86d74e61SPeter Brune 
249*86d74e61SPeter Brune .seealso: PetscLineSearchGetPreCheck(), PetscLineSearchSetPostCheck()
250*86d74e61SPeter Brune @*/
251*86d74e61SPeter Brune PetscErrorCode  PetscLineSearchGetPostCheck(PetscLineSearch linesearch, PetscLineSearchPostCheckFunc *func,void **ctx)
252*86d74e61SPeter Brune {
253*86d74e61SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
254*86d74e61SPeter Brune   if (func) *func = linesearch->ops->postcheckstep;
255*86d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
256*86d74e61SPeter Brune   PetscFunctionReturn(0);
257*86d74e61SPeter Brune }
258*86d74e61SPeter Brune 
259*86d74e61SPeter Brune 
260bf7f4e0aSPeter Brune #undef __FUNCT__
2616188f407SPeter Brune #define __FUNCT__ "PetscLineSearchPreCheck"
262f40b411bSPeter Brune /*@
2636188f407SPeter Brune    PetscLineSearchPreCheck - Prepares the line search for being applied.
264f40b411bSPeter Brune 
2656188f407SPeter Brune    Collective on PetscLineSearch
266f40b411bSPeter Brune 
267f40b411bSPeter Brune    Input Parameters:
268f40b411bSPeter Brune .  linesearch - The linesearch instance.
269f40b411bSPeter Brune 
270f40b411bSPeter Brune    Output Parameters:
271f40b411bSPeter Brune .  changed - Indicator if the pre-check has changed anything.
272f40b411bSPeter Brune 
273f40b411bSPeter Brune    Level: Beginner
274f40b411bSPeter Brune 
2756188f407SPeter Brune    .keywords: PetscLineSearch, Create
276f40b411bSPeter Brune 
2776188f407SPeter Brune    .seealso: PetscLineSearchPostCheck()
278f40b411bSPeter Brune @*/
2796188f407SPeter Brune PetscErrorCode PetscLineSearchPreCheck(PetscLineSearch linesearch, PetscBool * changed)
280bf7f4e0aSPeter Brune {
281bf7f4e0aSPeter Brune   PetscErrorCode ierr;
282bf7f4e0aSPeter Brune   PetscFunctionBegin;
283bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
284bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
285bf7f4e0aSPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed);CHKERRQ(ierr);
286bf7f4e0aSPeter Brune   }
287bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
288bf7f4e0aSPeter Brune }
289bf7f4e0aSPeter Brune 
290bf7f4e0aSPeter Brune #undef __FUNCT__
2916188f407SPeter Brune #define __FUNCT__ "PetscLineSearchPostCheck"
292f40b411bSPeter Brune /*@
2936188f407SPeter Brune    PetscLineSearchPostCheck - Prepares the line search for being applied.
294f40b411bSPeter Brune 
2956188f407SPeter Brune    Collective on PetscLineSearch
296f40b411bSPeter Brune 
297f40b411bSPeter Brune    Input Parameters:
298f40b411bSPeter Brune .  linesearch - The linesearch instance.
299f40b411bSPeter Brune 
300f40b411bSPeter Brune    Output Parameters:
301*86d74e61SPeter Brune +  changed_Y - Indicator if the solution has been changed.
302*86d74e61SPeter Brune -  changed_W - Indicator if the direction has been changed.
303f40b411bSPeter Brune 
304f40b411bSPeter Brune    Level: Intermediate
305f40b411bSPeter Brune 
3066188f407SPeter Brune    .keywords: PetscLineSearch, Create
307f40b411bSPeter Brune 
3086188f407SPeter Brune    .seealso: PetscLineSearchPreCheck()
309f40b411bSPeter Brune @*/
310*86d74e61SPeter Brune PetscErrorCode PetscLineSearchPostCheck(PetscLineSearch linesearch, PetscBool * changed_Y, PetscBool * changed_W)
311bf7f4e0aSPeter Brune {
312bf7f4e0aSPeter Brune   PetscErrorCode ierr;
313bf7f4e0aSPeter Brune   PetscFunctionBegin;
314bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
315bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
316bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
317*86d74e61SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, linesearch->vec_sol_new, changed_Y, changed_W);CHKERRQ(ierr);
318*86d74e61SPeter Brune   }
319*86d74e61SPeter Brune   PetscFunctionReturn(0);
320*86d74e61SPeter Brune }
321*86d74e61SPeter Brune 
322*86d74e61SPeter Brune 
323*86d74e61SPeter Brune #undef __FUNCT__
324*86d74e61SPeter Brune #define __FUNCT__ "PetscLineSearchPreCheckPicard"
325*86d74e61SPeter Brune /*@C
326*86d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
327*86d74e61SPeter Brune 
328*86d74e61SPeter Brune    Logically Collective
329*86d74e61SPeter Brune 
330*86d74e61SPeter Brune    Input Arguments:
331*86d74e61SPeter Brune +  linesearch - linesearch context
332*86d74e61SPeter Brune .  X - base state for this step
333*86d74e61SPeter Brune .  Y - initial correction
334*86d74e61SPeter Brune 
335*86d74e61SPeter Brune    Output Arguments:
336*86d74e61SPeter Brune +  Y - correction, possibly modified
337*86d74e61SPeter Brune -  changed - flag indicating that Y was modified
338*86d74e61SPeter Brune 
339*86d74e61SPeter Brune    Options Database Key:
340*86d74e61SPeter Brune +  -snes_ls_precheck_picard - activate this routine
341*86d74e61SPeter Brune -  -snes_ls_precheck_picard_angle - angle
342*86d74e61SPeter Brune 
343*86d74e61SPeter Brune    Level: advanced
344*86d74e61SPeter Brune 
345*86d74e61SPeter Brune    Notes:
346*86d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
347*86d74e61SPeter Brune 
348*86d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
349*86d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
350*86d74e61SPeter Brune    is generally not useful when using a Newton linearization.
351*86d74e61SPeter Brune 
352*86d74e61SPeter Brune    Reference:
353*86d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
354*86d74e61SPeter Brune 
355*86d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
356*86d74e61SPeter Brune @*/
357*86d74e61SPeter Brune PetscErrorCode PetscLineSearchPreCheckPicard(PetscLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
358*86d74e61SPeter Brune {
359*86d74e61SPeter Brune   PetscErrorCode ierr;
360*86d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
361*86d74e61SPeter Brune   Vec            Ylast;
362*86d74e61SPeter Brune   PetscScalar    dot;
363*86d74e61SPeter Brune   PetscInt       iter;
364*86d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
365*86d74e61SPeter Brune   SNES           snes;
366*86d74e61SPeter Brune 
367*86d74e61SPeter Brune   PetscFunctionBegin;
368*86d74e61SPeter Brune   ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
369*86d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
370*86d74e61SPeter Brune   if (!Ylast) {
371*86d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
372*86d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
373*86d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
374*86d74e61SPeter Brune   }
375*86d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
376*86d74e61SPeter Brune   if (iter < 2) {
377*86d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
378*86d74e61SPeter Brune     *changed = PETSC_FALSE;
379*86d74e61SPeter Brune     PetscFunctionReturn(0);
380*86d74e61SPeter Brune   }
381*86d74e61SPeter Brune 
382*86d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
383*86d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
384*86d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
385*86d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
386*86d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
387*86d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
388*86d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
389*86d74e61SPeter Brune     /* Modify the step Y */
390*86d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
391*86d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
392*86d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
393*86d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
394*86d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
395*86d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
396*86d74e61SPeter Brune     ierr = PetscInfo3(snes,"Angle %G degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,angle,alpha);CHKERRQ(ierr);
397*86d74e61SPeter Brune   } else {
398*86d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
399*86d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
400*86d74e61SPeter Brune     *changed = PETSC_FALSE;
401bf7f4e0aSPeter Brune   }
402bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
403bf7f4e0aSPeter Brune }
404bf7f4e0aSPeter Brune 
405bf7f4e0aSPeter Brune #undef __FUNCT__
4066188f407SPeter Brune #define __FUNCT__ "PetscLineSearchApply"
407f40b411bSPeter Brune /*@
4086188f407SPeter Brune    PetscLineSearchApply - Computes the line-search update
409f40b411bSPeter Brune 
4106188f407SPeter Brune    Collective on PetscLineSearch
411f40b411bSPeter Brune 
412f40b411bSPeter Brune    Input Parameters:
413f40b411bSPeter Brune +  linesearch - The linesearch instance.
414f40b411bSPeter Brune .  X - The current solution.
415f40b411bSPeter Brune .  F - The current function.
416f40b411bSPeter Brune .  fnorm - The current norm.
417f40b411bSPeter Brune .  Y - The search direction.
418f40b411bSPeter Brune 
419f40b411bSPeter Brune    Output Parameters:
420f40b411bSPeter Brune +  X - The new solution.
421f40b411bSPeter Brune .  F - The new function.
422f40b411bSPeter Brune -  fnorm - The new function norm.
423f40b411bSPeter Brune 
424f40b411bSPeter Brune    Level: Intermediate
425f40b411bSPeter Brune 
4266188f407SPeter Brune    .keywords: PetscLineSearch, Create
427f40b411bSPeter Brune 
4286188f407SPeter Brune    .seealso: PetscLineSearchCreate(), PetscLineSearchPreCheck(), PetscLineSearchPostCheck()
429f40b411bSPeter Brune @*/
4306188f407SPeter Brune PetscErrorCode PetscLineSearchApply(PetscLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
431bf7f4e0aSPeter Brune   PetscErrorCode ierr;
432bf7f4e0aSPeter Brune   PetscFunctionBegin;
433bf7f4e0aSPeter Brune 
434bf7f4e0aSPeter Brune   /* check the pointers */
4356188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
436bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
437bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
438bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
439bf7f4e0aSPeter Brune 
440bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
441bf7f4e0aSPeter Brune 
442bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
443bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
444bf7f4e0aSPeter Brune   linesearch->vec_func = F;
445bf7f4e0aSPeter Brune 
4466188f407SPeter Brune   ierr = PetscLineSearchSetUp(linesearch);CHKERRQ(ierr);
447bf7f4e0aSPeter Brune 
448bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
449bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
450bf7f4e0aSPeter Brune 
451bf7f4e0aSPeter Brune   if (fnorm) {
452bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
453bf7f4e0aSPeter Brune   } else {
454bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
455bf7f4e0aSPeter Brune   }
456bf7f4e0aSPeter Brune 
4576188f407SPeter Brune   ierr = PetscLogEventBegin(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
458bf7f4e0aSPeter Brune 
459bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
460bf7f4e0aSPeter Brune 
4616188f407SPeter Brune   ierr = PetscLogEventEnd(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
462bf7f4e0aSPeter Brune 
463bf7f4e0aSPeter Brune   if (fnorm)
464bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
465bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
466bf7f4e0aSPeter Brune }
467bf7f4e0aSPeter Brune 
468bf7f4e0aSPeter Brune #undef __FUNCT__
4696188f407SPeter Brune #define __FUNCT__ "PetscLineSearchDestroy"
470f40b411bSPeter Brune /*@
4716188f407SPeter Brune    PetscLineSearchDestroy - Destroys the line search instance.
472f40b411bSPeter Brune 
4736188f407SPeter Brune    Collective on PetscLineSearch
474f40b411bSPeter Brune 
475f40b411bSPeter Brune    Input Parameters:
476f40b411bSPeter Brune .  linesearch - The linesearch instance.
477f40b411bSPeter Brune 
478f40b411bSPeter Brune    Level: Intermediate
479f40b411bSPeter Brune 
4806188f407SPeter Brune    .keywords: PetscLineSearch, Create
481f40b411bSPeter Brune 
4826188f407SPeter Brune    .seealso: PetscLineSearchCreate(), PetscLineSearchReset()
483f40b411bSPeter Brune @*/
4846188f407SPeter Brune PetscErrorCode PetscLineSearchDestroy(PetscLineSearch * linesearch) {
485bf7f4e0aSPeter Brune   PetscErrorCode ierr;
486bf7f4e0aSPeter Brune   PetscFunctionBegin;
487bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
4886188f407SPeter Brune   PetscValidHeaderSpecific((*linesearch),PETSCLINESEARCH_CLASSID,1);
489bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
490bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
4916188f407SPeter Brune   ierr = PetscLineSearchReset(*linesearch);
492bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
493bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
494bf7f4e0aSPeter Brune   }
495bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
496e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
497bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
498bf7f4e0aSPeter Brune }
499bf7f4e0aSPeter Brune 
500bf7f4e0aSPeter Brune #undef __FUNCT__
5016188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetMonitor"
502f40b411bSPeter Brune /*@
5036188f407SPeter Brune    PetscLineSearchSetMonitor - Turns on/off printing useful things about the line search.
504bf7f4e0aSPeter Brune 
505bf7f4e0aSPeter Brune    Input Parameters:
506bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
507bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
508bf7f4e0aSPeter Brune 
509bf7f4e0aSPeter Brune    Logically Collective on SNES
510bf7f4e0aSPeter Brune 
511bf7f4e0aSPeter Brune    Options Database:
512f40b411bSPeter Brune .   -linesearch_monitor - enables the monitor.
513bf7f4e0aSPeter Brune 
514bf7f4e0aSPeter Brune    Level: intermediate
515bf7f4e0aSPeter Brune 
516bf7f4e0aSPeter Brune 
5176188f407SPeter Brune .seealso: PetscLineSearchGetMonitor()
518bf7f4e0aSPeter Brune @*/
5196188f407SPeter Brune PetscErrorCode  PetscLineSearchSetMonitor(PetscLineSearch linesearch, PetscBool flg)
520bf7f4e0aSPeter Brune {
521bf7f4e0aSPeter Brune 
522bf7f4e0aSPeter Brune   PetscErrorCode ierr;
523bf7f4e0aSPeter Brune   PetscFunctionBegin;
524bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
525bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
526bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
527bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
528bf7f4e0aSPeter Brune   }
529bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
530bf7f4e0aSPeter Brune }
531bf7f4e0aSPeter Brune 
532bf7f4e0aSPeter Brune #undef __FUNCT__
5336188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetMonitor"
534f40b411bSPeter Brune /*@
5356188f407SPeter Brune    PetscLineSearchGetMonitor - Gets the monitor instance for the line search
5366a388c36SPeter Brune 
537f40b411bSPeter Brune    Input Parameters:
538f40b411bSPeter Brune .  linesearch - linesearch context.
539f40b411bSPeter Brune 
540f40b411bSPeter Brune    Input Parameters:
541f40b411bSPeter Brune .  monitor - monitor context.
542f40b411bSPeter Brune 
543f40b411bSPeter Brune    Logically Collective on SNES
544f40b411bSPeter Brune 
545f40b411bSPeter Brune 
546f40b411bSPeter Brune    Options Database Keys:
547f40b411bSPeter Brune .   -linesearch_monitor - enables the monitor.
548f40b411bSPeter Brune 
549f40b411bSPeter Brune    Level: intermediate
550f40b411bSPeter Brune 
551f40b411bSPeter Brune 
5526188f407SPeter Brune .seealso: PetscLineSearchSetMonitor()
553f40b411bSPeter Brune @*/
5546188f407SPeter Brune PetscErrorCode  PetscLineSearchGetMonitor(PetscLineSearch linesearch, PetscViewer *monitor)
5556a388c36SPeter Brune {
5566a388c36SPeter Brune 
5576a388c36SPeter Brune   PetscFunctionBegin;
5586188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
5596a388c36SPeter Brune   if (monitor) {
5606a388c36SPeter Brune     PetscValidPointer(monitor, 2);
5616a388c36SPeter Brune     *monitor = linesearch->monitor;
5626a388c36SPeter Brune   }
5636a388c36SPeter Brune   PetscFunctionReturn(0);
5646a388c36SPeter Brune }
5656a388c36SPeter Brune 
5666a388c36SPeter Brune #undef __FUNCT__
5676188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetFromOptions"
568f40b411bSPeter Brune /*@
5696188f407SPeter Brune    PetscLineSearchSetFromOptions - Sets options for the line search
570f40b411bSPeter Brune 
571f40b411bSPeter Brune    Input Parameters:
572f40b411bSPeter Brune .  linesearch - linesearch context.
573f40b411bSPeter Brune 
574f40b411bSPeter Brune    Options Database Keys:
575f40b411bSPeter Brune + -linesearch_type - The Line search method
576f40b411bSPeter Brune . -linesearch_monitor - Print progress of line searches
577f40b411bSPeter Brune . -linesearch_damping - The linesearch damping parameter.
578f40b411bSPeter Brune . -linesearch_norms   - Turn on/off the linesearch norms
579f40b411bSPeter Brune . -linesearch_keeplambda - Keep the previous search length as the initial guess.
580f40b411bSPeter Brune - -linesearch_max_it - The number of iterations for iterative line searches.
581f40b411bSPeter Brune 
5826188f407SPeter Brune    Logically Collective on PetscLineSearch
583f40b411bSPeter Brune 
584f40b411bSPeter Brune    Level: intermediate
585f40b411bSPeter Brune 
586f40b411bSPeter Brune 
5876188f407SPeter Brune .seealso: PetscLineSearchCreate()
588f40b411bSPeter Brune @*/
5896188f407SPeter Brune PetscErrorCode PetscLineSearchSetFromOptions(PetscLineSearch linesearch) {
590bf7f4e0aSPeter Brune   PetscErrorCode ierr;
5916188f407SPeter Brune   const char     *deft = PETSCLINESEARCHBASIC;
592bf7f4e0aSPeter Brune   char           type[256];
593bf7f4e0aSPeter Brune   PetscBool      flg, set;
594bf7f4e0aSPeter Brune   PetscFunctionBegin;
5956188f407SPeter Brune   if (!PetscLineSearchRegisterAllCalled) {ierr = PetscLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
596bf7f4e0aSPeter Brune 
597bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
598bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
599bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
600bf7f4e0aSPeter Brune   }
6016188f407SPeter Brune   ierr = PetscOptionsList("-linesearch_type","Line-search method","PetscLineSearchSetType",PetscLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
602bf7f4e0aSPeter Brune   if (flg) {
6036188f407SPeter Brune     ierr = PetscLineSearchSetType(linesearch,type);CHKERRQ(ierr);
604bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
6056188f407SPeter Brune     ierr = PetscLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
606bf7f4e0aSPeter Brune   }
607bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
608bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
609bf7f4e0aSPeter Brune   }
610bf7f4e0aSPeter Brune 
6116188f407SPeter Brune   ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESPetscLineSearchSetMonitor",
612bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
6136188f407SPeter Brune   if (set) {ierr = PetscLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
614bf7f4e0aSPeter Brune 
6156188f407SPeter Brune   ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","PetscLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6166188f407SPeter Brune   ierr = PetscOptionsReal("-linesearch_rtol","Tolerance for iterative line search","PetscLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6176188f407SPeter Brune   ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","PetscLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
6186188f407SPeter Brune   ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","PetscLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
619bf7f4e0aSPeter Brune   ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
620bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
621bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
622bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
623bf7f4e0aSPeter Brune }
624bf7f4e0aSPeter Brune 
625bf7f4e0aSPeter Brune #undef __FUNCT__
6266188f407SPeter Brune #define __FUNCT__ "PetscLineSearchView"
627f40b411bSPeter Brune /*@
6286188f407SPeter Brune    PetscLineSearchView - Views useful information for the line search.
629f40b411bSPeter Brune 
630f40b411bSPeter Brune    Input Parameters:
631f40b411bSPeter Brune .  linesearch - linesearch context.
632f40b411bSPeter Brune 
6336188f407SPeter Brune    Logically Collective on PetscLineSearch
634f40b411bSPeter Brune 
635f40b411bSPeter Brune    Level: intermediate
636f40b411bSPeter Brune 
637f40b411bSPeter Brune 
6386188f407SPeter Brune .seealso: PetscLineSearchCreate()
639f40b411bSPeter Brune @*/
6406188f407SPeter Brune PetscErrorCode PetscLineSearchView(PetscLineSearch linesearch) {
641bf7f4e0aSPeter Brune   PetscFunctionBegin;
642f40b411bSPeter Brune 
643bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
644bf7f4e0aSPeter Brune }
645bf7f4e0aSPeter Brune 
646bf7f4e0aSPeter Brune #undef __FUNCT__
6476188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetType"
648f40b411bSPeter Brune /*@
6496188f407SPeter Brune    PetscLineSearchSetType - Sets the linesearch type
650f40b411bSPeter Brune 
651f40b411bSPeter Brune    Input Parameters:
652f40b411bSPeter Brune +  linesearch - linesearch context.
653f40b411bSPeter Brune -  type - The type of line search to be used
654f40b411bSPeter Brune 
6556188f407SPeter Brune    Logically Collective on PetscLineSearch
656f40b411bSPeter Brune 
657f40b411bSPeter Brune    Level: intermediate
658f40b411bSPeter Brune 
659f40b411bSPeter Brune 
6606188f407SPeter Brune .seealso: PetscLineSearchCreate()
661f40b411bSPeter Brune @*/
6626188f407SPeter Brune PetscErrorCode PetscLineSearchSetType(PetscLineSearch linesearch, const PetscLineSearchType type)
663bf7f4e0aSPeter Brune {
664bf7f4e0aSPeter Brune 
6656188f407SPeter Brune   PetscErrorCode ierr,(*r)(PetscLineSearch);
666bf7f4e0aSPeter Brune   PetscBool      match;
667bf7f4e0aSPeter Brune 
668bf7f4e0aSPeter Brune   PetscFunctionBegin;
6696188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
670bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
671bf7f4e0aSPeter Brune 
672bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
673bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
674bf7f4e0aSPeter Brune 
6756188f407SPeter Brune   ierr =  PetscFListFind(PetscLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
676bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
677bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
678bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
679bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
680bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
681bf7f4e0aSPeter Brune   }
6826188f407SPeter Brune   /* Reinitialize function pointers in PetscLineSearchOps structure */
683bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
684bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
685bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
686bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
687bf7f4e0aSPeter Brune 
688bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
689bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
690bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
691bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
692bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
693bf7f4e0aSPeter Brune   }
694bf7f4e0aSPeter Brune #endif
695bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
696bf7f4e0aSPeter Brune }
697bf7f4e0aSPeter Brune 
698bf7f4e0aSPeter Brune #undef __FUNCT__
6996188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetSNES"
700f40b411bSPeter Brune /*@
7016188f407SPeter Brune    PetscLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation
702f40b411bSPeter Brune 
703f40b411bSPeter Brune    Input Parameters:
704f40b411bSPeter Brune +  linesearch - linesearch context.
705f40b411bSPeter Brune -  snes - The snes instance
706f40b411bSPeter Brune 
707f40b411bSPeter Brune    Level: intermediate
708f40b411bSPeter Brune 
709f40b411bSPeter Brune 
7106188f407SPeter Brune .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
711f40b411bSPeter Brune @*/
7126188f407SPeter Brune PetscErrorCode  PetscLineSearchSetSNES(PetscLineSearch linesearch, SNES snes){
713bf7f4e0aSPeter Brune   PetscFunctionBegin;
7146188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
715bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
716bf7f4e0aSPeter Brune   linesearch->snes = snes;
717bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
718bf7f4e0aSPeter Brune }
719bf7f4e0aSPeter Brune 
720bf7f4e0aSPeter Brune #undef __FUNCT__
7216188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetSNES"
722f40b411bSPeter Brune /*@
7236188f407SPeter Brune    PetscLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation
724f40b411bSPeter Brune 
725f40b411bSPeter Brune    Input Parameters:
726f40b411bSPeter Brune .  linesearch - linesearch context.
727f40b411bSPeter Brune 
728f40b411bSPeter Brune    Output Parameters:
729f40b411bSPeter Brune .  snes - The snes instance
730f40b411bSPeter Brune 
731f40b411bSPeter Brune    Level: intermediate
732f40b411bSPeter Brune 
7336188f407SPeter Brune .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
734f40b411bSPeter Brune @*/
7356188f407SPeter Brune PetscErrorCode  PetscLineSearchGetSNES(PetscLineSearch linesearch, SNES *snes){
736bf7f4e0aSPeter Brune   PetscFunctionBegin;
7376188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
7386a388c36SPeter Brune   PetscValidPointer(snes, 2);
739bf7f4e0aSPeter Brune   *snes = linesearch->snes;
740bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
741bf7f4e0aSPeter Brune }
742bf7f4e0aSPeter Brune 
7436a388c36SPeter Brune #undef __FUNCT__
7446188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetLambda"
745f40b411bSPeter Brune /*@
7466188f407SPeter Brune    PetscLineSearchGetLambda - Gets the last linesearch steplength discovered.
747f40b411bSPeter Brune 
748f40b411bSPeter Brune    Input Parameters:
749f40b411bSPeter Brune .  linesearch - linesearch context.
750f40b411bSPeter Brune 
751f40b411bSPeter Brune    Output Parameters:
752f40b411bSPeter Brune .  lambda - The last steplength.
753f40b411bSPeter Brune 
754f40b411bSPeter Brune    Level: intermediate
755f40b411bSPeter Brune 
7566188f407SPeter Brune .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
757f40b411bSPeter Brune @*/
7586188f407SPeter Brune PetscErrorCode  PetscLineSearchGetLambda(PetscLineSearch linesearch,PetscReal *lambda)
7596a388c36SPeter Brune {
7606a388c36SPeter Brune   PetscFunctionBegin;
7616188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
7626a388c36SPeter Brune   PetscValidPointer(lambda, 2);
7636a388c36SPeter Brune   *lambda = linesearch->lambda;
7646a388c36SPeter Brune   PetscFunctionReturn(0);
7656a388c36SPeter Brune }
7666a388c36SPeter Brune 
7676a388c36SPeter Brune #undef __FUNCT__
7686188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetLambda"
769f40b411bSPeter Brune /*@
7706188f407SPeter Brune    PetscLineSearchSetLambda - Sets the linesearch steplength.
771f40b411bSPeter Brune 
772f40b411bSPeter Brune    Input Parameters:
773f40b411bSPeter Brune +  linesearch - linesearch context.
774f40b411bSPeter Brune -  lambda - The last steplength.
775f40b411bSPeter Brune 
776f40b411bSPeter Brune    Level: intermediate
777f40b411bSPeter Brune 
7786188f407SPeter Brune .seealso: PetscLineSearchGetLambda()
779f40b411bSPeter Brune @*/
7806188f407SPeter Brune PetscErrorCode  PetscLineSearchSetLambda(PetscLineSearch linesearch, PetscReal lambda)
7816a388c36SPeter Brune {
7826a388c36SPeter Brune   PetscFunctionBegin;
7836188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
7846a388c36SPeter Brune   linesearch->lambda = lambda;
7856a388c36SPeter Brune   PetscFunctionReturn(0);
7866a388c36SPeter Brune }
7876a388c36SPeter Brune 
7886a388c36SPeter Brune #undef  __FUNCT__
7896188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetTolerances"
790f40b411bSPeter Brune /*@
7916188f407SPeter Brune    PetscLineSearchGetTolerances - Gets the tolerances for the method
792f40b411bSPeter Brune 
793f40b411bSPeter Brune    Input Parameters:
794f40b411bSPeter Brune .  linesearch - linesearch context.
795f40b411bSPeter Brune 
796f40b411bSPeter Brune    Output Parameters:
797516fe3c3SPeter Brune +  steptol - The minimum steplength
798516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
799516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
800516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
801516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
802f40b411bSPeter Brune 
803f40b411bSPeter Brune 
804516fe3c3SPeter Brune    Level: advanced
805516fe3c3SPeter Brune 
8066188f407SPeter Brune .seealso: PetscLineSearchSetTolerances()
807f40b411bSPeter Brune @*/
8086188f407SPeter Brune PetscErrorCode  PetscLineSearchGetTolerances(PetscLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
8096a388c36SPeter Brune {
8106a388c36SPeter Brune   PetscFunctionBegin;
8116188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
812516fe3c3SPeter Brune   if (steptol) {
8136a388c36SPeter Brune     PetscValidPointer(steptol, 2);
8146a388c36SPeter Brune     *steptol = linesearch->steptol;
815516fe3c3SPeter Brune   }
816516fe3c3SPeter Brune   if (maxstep) {
817516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
818516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
819516fe3c3SPeter Brune   }
820516fe3c3SPeter Brune   if (rtol) {
821516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
822516fe3c3SPeter Brune     *rtol = linesearch->rtol;
823516fe3c3SPeter Brune   }
824516fe3c3SPeter Brune   if (atol) {
825516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
826516fe3c3SPeter Brune     *atol = linesearch->atol;
827516fe3c3SPeter Brune   }
828516fe3c3SPeter Brune   if (ltol) {
829516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
830516fe3c3SPeter Brune     *ltol = linesearch->ltol;
831516fe3c3SPeter Brune   }
832516fe3c3SPeter Brune   if (max_its) {
833516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
834516fe3c3SPeter Brune     *max_its = linesearch->max_its;
835516fe3c3SPeter Brune   }
8366a388c36SPeter Brune   PetscFunctionReturn(0);
8376a388c36SPeter Brune }
8386a388c36SPeter Brune 
8396a388c36SPeter Brune #undef  __FUNCT__
8406188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetTolerances"
841f40b411bSPeter Brune /*@
8426188f407SPeter Brune    PetscLineSearchSetTolerances - Sets the tolerances for the method
843f40b411bSPeter Brune 
844f40b411bSPeter Brune    Input Parameters:
845516fe3c3SPeter Brune +  linesearch - linesearch context.
846516fe3c3SPeter Brune .  steptol - The minimum steplength
847516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
848516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
849516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
850516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
851f40b411bSPeter Brune 
852f40b411bSPeter Brune 
853516fe3c3SPeter Brune    Level: advanced
854516fe3c3SPeter Brune 
8556188f407SPeter Brune .seealso: PetscLineSearchGetTolerances()
856f40b411bSPeter Brune @*/
8576188f407SPeter Brune PetscErrorCode  PetscLineSearchSetTolerances(PetscLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
8586a388c36SPeter Brune {
8596a388c36SPeter Brune   PetscFunctionBegin;
8606188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
8616a388c36SPeter Brune   linesearch->steptol = steptol;
862516fe3c3SPeter Brune   linesearch->maxstep = maxstep;
863516fe3c3SPeter Brune   linesearch->rtol = rtol;
864516fe3c3SPeter Brune   linesearch->atol = atol;
865516fe3c3SPeter Brune   linesearch->ltol = ltol;
866516fe3c3SPeter Brune   linesearch->max_its = max_its;
8676a388c36SPeter Brune   PetscFunctionReturn(0);
8686a388c36SPeter Brune }
8696a388c36SPeter Brune 
870516fe3c3SPeter Brune 
8716a388c36SPeter Brune #undef __FUNCT__
8726188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetDamping"
873f40b411bSPeter Brune /*@
8746188f407SPeter Brune    PetscLineSearchGetDamping - Gets the line search damping parameter.
875f40b411bSPeter Brune 
876f40b411bSPeter Brune    Input Parameters:
877f40b411bSPeter Brune .  linesearch - linesearch context.
878f40b411bSPeter Brune 
879f40b411bSPeter Brune    Output Parameters:
880f40b411bSPeter Brune .  damping - The damping parameter.
881f40b411bSPeter Brune 
882f40b411bSPeter Brune    Level: intermediate
883f40b411bSPeter Brune 
8846188f407SPeter Brune .seealso: PetscLineSearchGetStepTolerance()
885f40b411bSPeter Brune @*/
886f40b411bSPeter Brune 
8876188f407SPeter Brune PetscErrorCode  PetscLineSearchGetDamping(PetscLineSearch linesearch,PetscReal *damping)
8886a388c36SPeter Brune {
8896a388c36SPeter Brune   PetscFunctionBegin;
8906188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
8916a388c36SPeter Brune   PetscValidPointer(damping, 2);
8926a388c36SPeter Brune   *damping = linesearch->damping;
8936a388c36SPeter Brune   PetscFunctionReturn(0);
8946a388c36SPeter Brune }
8956a388c36SPeter Brune 
8966a388c36SPeter Brune #undef __FUNCT__
8976188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetDamping"
898f40b411bSPeter Brune /*@
8996188f407SPeter Brune    PetscLineSearchSetDamping - Sets the line search damping paramter.
900f40b411bSPeter Brune 
901f40b411bSPeter Brune    Input Parameters:
902f40b411bSPeter Brune .  linesearch - linesearch context.
903f40b411bSPeter Brune .  damping - The damping parameter.
904f40b411bSPeter Brune 
905f40b411bSPeter Brune    Level: intermediate
906f40b411bSPeter Brune 
9076188f407SPeter Brune .seealso: PetscLineSearchGetDamping()
908f40b411bSPeter Brune @*/
9096188f407SPeter Brune PetscErrorCode  PetscLineSearchSetDamping(PetscLineSearch linesearch,PetscReal damping)
9106a388c36SPeter Brune {
9116a388c36SPeter Brune   PetscFunctionBegin;
9126188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
9136a388c36SPeter Brune   linesearch->damping = damping;
9146a388c36SPeter Brune   PetscFunctionReturn(0);
9156a388c36SPeter Brune }
9166a388c36SPeter Brune 
9176a388c36SPeter Brune #undef __FUNCT__
9186188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetNorms"
919f40b411bSPeter Brune /*@
9206188f407SPeter Brune    PetscLineSearchGetNorms - Gets the norms for for X, Y, and F.
921f40b411bSPeter Brune 
922f40b411bSPeter Brune    Input Parameters:
923f40b411bSPeter Brune .  linesearch - linesearch context.
924f40b411bSPeter Brune 
925f40b411bSPeter Brune    Output Parameters:
926f40b411bSPeter Brune +  xnorm - The norm of the current solution
927f40b411bSPeter Brune .  fnorm - The norm of the current function
928f40b411bSPeter Brune -  ynorm - The norm of the current update
929f40b411bSPeter Brune 
930f40b411bSPeter Brune    Level: intermediate
931f40b411bSPeter Brune 
9326188f407SPeter Brune .seealso: PetscLineSearchSetNorms() PetscLineSearchGetVecs()
933f40b411bSPeter Brune @*/
9346188f407SPeter Brune PetscErrorCode  PetscLineSearchGetNorms(PetscLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
935bf7f4e0aSPeter Brune {
936bf7f4e0aSPeter Brune   PetscFunctionBegin;
9376188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
938bf7f4e0aSPeter Brune   if (xnorm) {
939bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
940bf7f4e0aSPeter Brune   }
941bf7f4e0aSPeter Brune   if (fnorm) {
942bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
943bf7f4e0aSPeter Brune   }
944bf7f4e0aSPeter Brune   if (ynorm) {
945bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
946bf7f4e0aSPeter Brune   }
947bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
948bf7f4e0aSPeter Brune }
949bf7f4e0aSPeter Brune 
950e7058c64SPeter Brune #undef __FUNCT__
9516188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetNorms"
952f40b411bSPeter Brune /*@
9536188f407SPeter Brune    PetscLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
954f40b411bSPeter Brune 
955f40b411bSPeter Brune    Input Parameters:
956f40b411bSPeter Brune +  linesearch - linesearch context.
957f40b411bSPeter Brune .  xnorm - The norm of the current solution
958f40b411bSPeter Brune .  fnorm - The norm of the current function
959f40b411bSPeter Brune -  ynorm - The norm of the current update
960f40b411bSPeter Brune 
961f40b411bSPeter Brune    Level: intermediate
962f40b411bSPeter Brune 
9636188f407SPeter Brune .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs()
964f40b411bSPeter Brune @*/
9656188f407SPeter Brune PetscErrorCode  PetscLineSearchSetNorms(PetscLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
9666a388c36SPeter Brune {
9676a388c36SPeter Brune   PetscFunctionBegin;
9686188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
9696a388c36SPeter Brune   linesearch->xnorm = xnorm;
9706a388c36SPeter Brune   linesearch->fnorm = fnorm;
9716a388c36SPeter Brune   linesearch->ynorm = ynorm;
9726a388c36SPeter Brune   PetscFunctionReturn(0);
9736a388c36SPeter Brune }
9746a388c36SPeter Brune 
9756a388c36SPeter Brune #undef __FUNCT__
9766188f407SPeter Brune #define __FUNCT__ "PetscLineSearchComputeNorms"
977f40b411bSPeter Brune /*@
9786188f407SPeter Brune    PetscLineSearchComputeNorms - Computes the norms of X, F, and Y.
979f40b411bSPeter Brune 
980f40b411bSPeter Brune    Input Parameters:
981f40b411bSPeter Brune .  linesearch - linesearch context.
982f40b411bSPeter Brune 
983f40b411bSPeter Brune    Options Database Keys:
984f40b411bSPeter Brune .   -linesearch_norms - turn norm computation on or off.
985f40b411bSPeter Brune 
986f40b411bSPeter Brune    Level: intermediate
987f40b411bSPeter Brune 
9886188f407SPeter Brune .seealso: PetscLineSearchGetNorms, PetscLineSearchSetNorms()
989f40b411bSPeter Brune @*/
9906188f407SPeter Brune PetscErrorCode PetscLineSearchComputeNorms(PetscLineSearch linesearch)
9916a388c36SPeter Brune {
9926a388c36SPeter Brune   PetscErrorCode ierr;
9936a388c36SPeter Brune   PetscFunctionBegin;
9946a388c36SPeter Brune   if (linesearch->norms) {
9956a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
9966a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
9976a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
9986a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
9996a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
10006a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
10016a388c36SPeter Brune   }
10026a388c36SPeter Brune   PetscFunctionReturn(0);
10036a388c36SPeter Brune }
10046a388c36SPeter Brune 
10056a388c36SPeter Brune #undef __FUNCT__
10066188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetVecs"
1007f40b411bSPeter Brune /*@
10086188f407SPeter Brune    PetscLineSearchGetVecs - Gets the vectors from the PetscLineSearch context
1009f40b411bSPeter Brune 
1010f40b411bSPeter Brune    Input Parameters:
1011f40b411bSPeter Brune .  linesearch - linesearch context.
1012f40b411bSPeter Brune 
1013f40b411bSPeter Brune    Output Parameters:
1014f40b411bSPeter Brune +  X - The old solution
1015f40b411bSPeter Brune .  F - The old function
1016f40b411bSPeter Brune .  Y - The search direction
1017f40b411bSPeter Brune .  W - The new solution
1018f40b411bSPeter Brune -  G - The new function
1019f40b411bSPeter Brune 
1020f40b411bSPeter Brune    Level: intermediate
1021f40b411bSPeter Brune 
10226188f407SPeter Brune .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs()
1023f40b411bSPeter Brune @*/
10246188f407SPeter Brune PetscErrorCode PetscLineSearchGetVecs(PetscLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
10256a388c36SPeter Brune   PetscFunctionBegin;
10266188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
10276a388c36SPeter Brune   if (X) {
10286a388c36SPeter Brune     PetscValidPointer(X, 2);
10296a388c36SPeter Brune     *X = linesearch->vec_sol;
10306a388c36SPeter Brune   }
10316a388c36SPeter Brune   if (F) {
10326a388c36SPeter Brune     PetscValidPointer(F, 3);
10336a388c36SPeter Brune     *F = linesearch->vec_func;
10346a388c36SPeter Brune   }
10356a388c36SPeter Brune   if (Y) {
10366a388c36SPeter Brune     PetscValidPointer(Y, 4);
10376a388c36SPeter Brune     *Y = linesearch->vec_update;
10386a388c36SPeter Brune   }
10396a388c36SPeter Brune   if (W) {
10406a388c36SPeter Brune     PetscValidPointer(W, 5);
10416a388c36SPeter Brune     *W = linesearch->vec_sol_new;
10426a388c36SPeter Brune   }
10436a388c36SPeter Brune   if (G) {
10446a388c36SPeter Brune     PetscValidPointer(G, 6);
10456a388c36SPeter Brune     *G = linesearch->vec_func_new;
10466a388c36SPeter Brune   }
10476a388c36SPeter Brune 
10486a388c36SPeter Brune   PetscFunctionReturn(0);
10496a388c36SPeter Brune }
10506a388c36SPeter Brune 
10516a388c36SPeter Brune #undef __FUNCT__
10526188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetVecs"
1053f40b411bSPeter Brune /*@
10546188f407SPeter Brune    PetscLineSearchSetVecs - Sets the vectors on the PetscLineSearch context
1055f40b411bSPeter Brune 
1056f40b411bSPeter Brune    Input Parameters:
1057f40b411bSPeter Brune +  linesearch - linesearch context.
1058f40b411bSPeter Brune .  X - The old solution
1059f40b411bSPeter Brune .  F - The old function
1060f40b411bSPeter Brune .  Y - The search direction
1061f40b411bSPeter Brune .  W - The new solution
1062f40b411bSPeter Brune -  G - The new function
1063f40b411bSPeter Brune 
1064f40b411bSPeter Brune    Level: intermediate
1065f40b411bSPeter Brune 
10666188f407SPeter Brune .seealso: PetscLineSearchSetNorms(), PetscLineSearchGetVecs()
1067f40b411bSPeter Brune @*/
10686188f407SPeter Brune PetscErrorCode PetscLineSearchSetVecs(PetscLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
10696a388c36SPeter Brune   PetscFunctionBegin;
10706188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
10716a388c36SPeter Brune   if (X) {
10726a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
10736a388c36SPeter Brune     linesearch->vec_sol = X;
10746a388c36SPeter Brune   }
10756a388c36SPeter Brune   if (F) {
10766a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
10776a388c36SPeter Brune     linesearch->vec_func = F;
10786a388c36SPeter Brune   }
10796a388c36SPeter Brune   if (Y) {
10806a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
10816a388c36SPeter Brune     linesearch->vec_update = Y;
10826a388c36SPeter Brune   }
10836a388c36SPeter Brune   if (W) {
10846a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
10856a388c36SPeter Brune     linesearch->vec_sol_new = W;
10866a388c36SPeter Brune   }
10876a388c36SPeter Brune   if (G) {
10886a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
10896a388c36SPeter Brune     linesearch->vec_func_new = G;
10906a388c36SPeter Brune   }
10916a388c36SPeter Brune 
10926a388c36SPeter Brune   PetscFunctionReturn(0);
10936a388c36SPeter Brune }
10946a388c36SPeter Brune 
10956a388c36SPeter Brune #undef __FUNCT__
10966188f407SPeter Brune #define __FUNCT__ "PetscLineSearchAppendOptionsPrefix"
1097e7058c64SPeter Brune /*@C
10986188f407SPeter Brune    PetscLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1099e7058c64SPeter Brune    SNES options in the database.
1100e7058c64SPeter Brune 
1101e7058c64SPeter Brune    Logically Collective on SNES
1102e7058c64SPeter Brune 
1103e7058c64SPeter Brune    Input Parameters:
1104e7058c64SPeter Brune +  snes - the SNES context
1105e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1106e7058c64SPeter Brune 
1107e7058c64SPeter Brune    Notes:
1108e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1109e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1110e7058c64SPeter Brune 
1111e7058c64SPeter Brune    Level: advanced
1112e7058c64SPeter Brune 
11136188f407SPeter Brune .keywords: PetscLineSearch, append, options, prefix, database
1114e7058c64SPeter Brune 
1115e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1116e7058c64SPeter Brune @*/
11176188f407SPeter Brune PetscErrorCode  PetscLineSearchAppendOptionsPrefix(PetscLineSearch linesearch,const char prefix[])
1118e7058c64SPeter Brune {
1119e7058c64SPeter Brune   PetscErrorCode ierr;
1120e7058c64SPeter Brune 
1121e7058c64SPeter Brune   PetscFunctionBegin;
11226188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
1123e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1124e7058c64SPeter Brune   PetscFunctionReturn(0);
1125e7058c64SPeter Brune }
1126e7058c64SPeter Brune 
1127e7058c64SPeter Brune #undef __FUNCT__
11286188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetOptionsPrefix"
1129e7058c64SPeter Brune /*@C
11306188f407SPeter Brune    PetscLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
11316188f407SPeter Brune    PetscLineSearch options in the database.
1132e7058c64SPeter Brune 
1133e7058c64SPeter Brune    Not Collective
1134e7058c64SPeter Brune 
1135e7058c64SPeter Brune    Input Parameter:
1136e7058c64SPeter Brune .  snes - the SNES context
1137e7058c64SPeter Brune 
1138e7058c64SPeter Brune    Output Parameter:
1139e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1140e7058c64SPeter Brune 
1141e7058c64SPeter Brune    Notes: On the fortran side, the user should pass in a string 'prefix' of
1142e7058c64SPeter Brune    sufficient length to hold the prefix.
1143e7058c64SPeter Brune 
1144e7058c64SPeter Brune    Level: advanced
1145e7058c64SPeter Brune 
11466188f407SPeter Brune .keywords: PetscLineSearch, get, options, prefix, database
1147e7058c64SPeter Brune 
1148e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1149e7058c64SPeter Brune @*/
11506188f407SPeter Brune PetscErrorCode  PetscLineSearchGetOptionsPrefix(PetscLineSearch linesearch,const char *prefix[])
1151e7058c64SPeter Brune {
1152e7058c64SPeter Brune   PetscErrorCode ierr;
1153e7058c64SPeter Brune 
1154e7058c64SPeter Brune   PetscFunctionBegin;
11556188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
1156e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1157e7058c64SPeter Brune   PetscFunctionReturn(0);
1158e7058c64SPeter Brune }
1159bf7f4e0aSPeter Brune 
1160bf7f4e0aSPeter Brune #undef __FUNCT__
11616188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetWork"
1162f40b411bSPeter Brune /*@
11636188f407SPeter Brune    PetscLineSearchGetWork - Gets work vectors for the line search.
1164f40b411bSPeter Brune 
1165f40b411bSPeter Brune    Input Parameter:
11666188f407SPeter Brune +  linesearch - the PetscLineSearch context
1167f40b411bSPeter Brune -  nwork - the number of work vectors
1168f40b411bSPeter Brune 
1169f40b411bSPeter Brune    Level: developer
1170f40b411bSPeter Brune 
11716188f407SPeter Brune .keywords: PetscLineSearch, work, vector
1172f40b411bSPeter Brune 
1173f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1174f40b411bSPeter Brune @*/
11756188f407SPeter Brune PetscErrorCode  PetscLineSearchGetWork(PetscLineSearch linesearch, PetscInt nwork)
1176bf7f4e0aSPeter Brune {
1177bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1178bf7f4e0aSPeter Brune   PetscFunctionBegin;
1179bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1180bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1181bf7f4e0aSPeter Brune   } else {
1182bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1183bf7f4e0aSPeter Brune   }
1184bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1185bf7f4e0aSPeter Brune }
1186bf7f4e0aSPeter Brune 
11876a388c36SPeter Brune 
1188bf7f4e0aSPeter Brune #undef __FUNCT__
11896188f407SPeter Brune #define __FUNCT__ "PetscLineSearchGetSuccess"
1190f40b411bSPeter Brune /*@
11916188f407SPeter Brune    PetscLineSearchGetSuccess - Gets the success/failure status of the last line search application
1192f40b411bSPeter Brune 
1193f40b411bSPeter Brune    Input Parameters:
1194f40b411bSPeter Brune .  linesearch - linesearch context.
1195f40b411bSPeter Brune 
1196f40b411bSPeter Brune    Output Parameters:
1197f40b411bSPeter Brune .  success - The success or failure status.
1198f40b411bSPeter Brune 
1199f40b411bSPeter Brune    Level: intermediate
1200f40b411bSPeter Brune 
12016188f407SPeter Brune .seealso: PetscLineSearchSetSuccess()
1202f40b411bSPeter Brune @*/
12036188f407SPeter Brune PetscErrorCode  PetscLineSearchGetSuccess(PetscLineSearch linesearch, PetscBool *success)
1204bf7f4e0aSPeter Brune {
1205bf7f4e0aSPeter Brune   PetscFunctionBegin;
12066188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
12076a388c36SPeter Brune   PetscValidPointer(success, 2);
1208bf7f4e0aSPeter Brune   if (success) {
1209bf7f4e0aSPeter Brune     *success = linesearch->success;
1210bf7f4e0aSPeter Brune   }
1211bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1212bf7f4e0aSPeter Brune }
1213bf7f4e0aSPeter Brune 
1214bf7f4e0aSPeter Brune #undef __FUNCT__
12156188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetSuccess"
1216f40b411bSPeter Brune /*@
12176188f407SPeter Brune    PetscLineSearchSetSuccess - Sets the success/failure status of the last line search application
1218f40b411bSPeter Brune 
1219f40b411bSPeter Brune    Input Parameters:
1220f40b411bSPeter Brune +  linesearch - linesearch context.
1221f40b411bSPeter Brune -  success - The success or failure status.
1222f40b411bSPeter Brune 
1223f40b411bSPeter Brune    Level: intermediate
1224f40b411bSPeter Brune 
12256188f407SPeter Brune .seealso: PetscLineSearchGetSuccess()
1226f40b411bSPeter Brune @*/
12276188f407SPeter Brune PetscErrorCode  PetscLineSearchSetSuccess(PetscLineSearch linesearch, PetscBool success)
12286a388c36SPeter Brune {
12296188f407SPeter Brune   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
12306a388c36SPeter Brune   PetscFunctionBegin;
12316a388c36SPeter Brune   linesearch->success = success;
12326a388c36SPeter Brune   PetscFunctionReturn(0);
12336a388c36SPeter Brune }
12346a388c36SPeter Brune 
12356a388c36SPeter Brune #undef __FUNCT__
12366188f407SPeter Brune #define __FUNCT__ "PetscLineSearchRegister"
1237bf7f4e0aSPeter Brune /*@C
12386188f407SPeter Brune   PetscLineSearchRegister - See PetscLineSearchRegisterDynamic()
1239bf7f4e0aSPeter Brune 
1240bf7f4e0aSPeter Brune   Level: advanced
1241bf7f4e0aSPeter Brune @*/
12426188f407SPeter Brune PetscErrorCode  PetscLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PetscLineSearch))
1243bf7f4e0aSPeter Brune {
1244bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1245bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1246bf7f4e0aSPeter Brune 
1247bf7f4e0aSPeter Brune   PetscFunctionBegin;
1248bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
12496188f407SPeter Brune   ierr = PetscFListAdd(&PetscLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1250bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1251bf7f4e0aSPeter Brune }
1252