xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 6a388c361ff2893f4363aecae9e02aa0b4e91a18)
1bf7f4e0aSPeter Brune #include <private/linesearchimpl.h> /*I "petsclinesearch.h" I*/
2bf7f4e0aSPeter Brune 
3bf7f4e0aSPeter Brune PetscBool  LineSearchRegisterAllCalled = PETSC_FALSE;
4bf7f4e0aSPeter Brune PetscFList LineSearchList              = PETSC_NULL;
5bf7f4e0aSPeter Brune 
6bf7f4e0aSPeter Brune PetscClassId   LineSearch_CLASSID;
7bf7f4e0aSPeter Brune PetscLogEvent  LineSearch_Apply;
8bf7f4e0aSPeter Brune 
9bf7f4e0aSPeter Brune #undef __FUNCT__
10bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchCreate"
11bf7f4e0aSPeter Brune PetscErrorCode LineSearchCreate(MPI_Comm comm, LineSearch * outlinesearch) {
12bf7f4e0aSPeter Brune   PetscErrorCode ierr;
13bf7f4e0aSPeter Brune   LineSearch     linesearch;
14bf7f4e0aSPeter Brune   PetscFunctionBegin;
15bf7f4e0aSPeter Brune   ierr = PetscHeaderCreate(linesearch, _p_LineSearch,struct _LineSearchOps,LineSearch_CLASSID, 0,
16bf7f4e0aSPeter Brune                            "LineSearch","Line-search method","LineSearch",comm,LineSearchDestroy,LineSearchView);CHKERRQ(ierr);
17bf7f4e0aSPeter Brune 
18bf7f4e0aSPeter Brune   linesearch->ops->precheckstep = PETSC_NULL;
19bf7f4e0aSPeter Brune   linesearch->ops->postcheckstep = PETSC_NULL;
20bf7f4e0aSPeter Brune 
21bf7f4e0aSPeter Brune   linesearch->lambda        = 1.0;
22bf7f4e0aSPeter Brune   linesearch->fnorm         = 1.0;
23bf7f4e0aSPeter Brune   linesearch->ynorm         = 1.0;
24bf7f4e0aSPeter Brune   linesearch->xnorm         = 1.0;
25bf7f4e0aSPeter Brune   linesearch->success       = PETSC_TRUE;
26bf7f4e0aSPeter Brune   linesearch->norms         = PETSC_TRUE;
27bf7f4e0aSPeter Brune   linesearch->keeplambda    = PETSC_FALSE;
28bf7f4e0aSPeter Brune   linesearch->damping       = 1.0;
29bf7f4e0aSPeter Brune   linesearch->maxstep       = 1e8;
30bf7f4e0aSPeter Brune   linesearch->steptol       = 1e-12;
31bf7f4e0aSPeter Brune   linesearch->precheckctx   = PETSC_NULL;
32bf7f4e0aSPeter Brune   linesearch->postcheckctx  = PETSC_NULL;
33bf7f4e0aSPeter Brune   linesearch->max_its       = 1;
34bf7f4e0aSPeter Brune   linesearch->setupcalled   = PETSC_FALSE;
35bf7f4e0aSPeter Brune   *outlinesearch            = linesearch;
36bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
37bf7f4e0aSPeter Brune }
38bf7f4e0aSPeter Brune 
39bf7f4e0aSPeter Brune #undef __FUNCT__
40bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetUp"
41bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetUp(LineSearch linesearch) {
42bf7f4e0aSPeter Brune   PetscErrorCode ierr;
43bf7f4e0aSPeter Brune   PetscFunctionBegin;
44bf7f4e0aSPeter Brune 
45bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
46bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,LINESEARCHBASIC);CHKERRQ(ierr);
47bf7f4e0aSPeter Brune   }
48bf7f4e0aSPeter Brune 
49bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
50bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
51bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_func, &linesearch->vec_func_new);CHKERRQ(ierr);
52bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
53bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
54bf7f4e0aSPeter Brune     }
55bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
56bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
57bf7f4e0aSPeter Brune   }
58bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
59bf7f4e0aSPeter Brune }
60bf7f4e0aSPeter Brune 
61bf7f4e0aSPeter Brune #undef __FUNCT__
62bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchReset"
63bf7f4e0aSPeter Brune PetscErrorCode LineSearchReset(LineSearch linesearch) {
64bf7f4e0aSPeter Brune   PetscErrorCode ierr;
65bf7f4e0aSPeter Brune   PetscFunctionBegin;
66bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
67bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
68bf7f4e0aSPeter Brune   }
69bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
70bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
71bf7f4e0aSPeter Brune 
72bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
73bf7f4e0aSPeter Brune   linesearch->nwork = 0;
74bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
75bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
76bf7f4e0aSPeter Brune }
77bf7f4e0aSPeter Brune 
78bf7f4e0aSPeter Brune #undef __FUNCT__
79bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchPreCheck"
80bf7f4e0aSPeter Brune PetscErrorCode LineSearchPreCheck(LineSearch linesearch, PetscBool * changed)
81bf7f4e0aSPeter Brune {
82bf7f4e0aSPeter Brune   PetscErrorCode ierr;
83bf7f4e0aSPeter Brune   PetscFunctionBegin;
84bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
85bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
86bf7f4e0aSPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed);CHKERRQ(ierr);
87bf7f4e0aSPeter Brune   }
88bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
89bf7f4e0aSPeter Brune }
90bf7f4e0aSPeter Brune 
91bf7f4e0aSPeter Brune #undef __FUNCT__
92bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchPostCheck"
93bf7f4e0aSPeter Brune PetscErrorCode LineSearchPostCheck(LineSearch linesearch, PetscBool * changed_W, PetscBool * changed_Y)
94bf7f4e0aSPeter Brune {
95bf7f4e0aSPeter Brune   PetscErrorCode ierr;
96bf7f4e0aSPeter Brune   PetscFunctionBegin;
97bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
98bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
99bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
100bf7f4e0aSPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_sol_new, linesearch->vec_update, changed_W, changed_Y);CHKERRQ(ierr);
101bf7f4e0aSPeter Brune   }
102bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
103bf7f4e0aSPeter Brune }
104bf7f4e0aSPeter Brune 
105bf7f4e0aSPeter Brune #undef __FUNCT__
106bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchApply"
107bf7f4e0aSPeter Brune PetscErrorCode LineSearchApply(LineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
108bf7f4e0aSPeter Brune   PetscErrorCode ierr;
109bf7f4e0aSPeter Brune   PetscFunctionBegin;
110bf7f4e0aSPeter Brune 
111bf7f4e0aSPeter Brune   /* check the pointers */
112bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
113bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
114bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
115bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
116bf7f4e0aSPeter Brune 
117bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
118bf7f4e0aSPeter Brune 
119bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
120bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
121bf7f4e0aSPeter Brune   linesearch->vec_func = F;
122bf7f4e0aSPeter Brune 
123bf7f4e0aSPeter Brune   ierr = LineSearchSetUp(linesearch);CHKERRQ(ierr);
124bf7f4e0aSPeter Brune 
125bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
126bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
127bf7f4e0aSPeter Brune 
128bf7f4e0aSPeter Brune   if (fnorm) {
129bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
130bf7f4e0aSPeter Brune   } else {
131bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
132bf7f4e0aSPeter Brune   }
133bf7f4e0aSPeter Brune 
134bf7f4e0aSPeter Brune   ierr = PetscLogEventBegin(LineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
135bf7f4e0aSPeter Brune 
136bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
137bf7f4e0aSPeter Brune 
138bf7f4e0aSPeter Brune   ierr = PetscLogEventEnd(LineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
139bf7f4e0aSPeter Brune 
140bf7f4e0aSPeter Brune   if (fnorm)
141bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
142bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
143bf7f4e0aSPeter Brune }
144bf7f4e0aSPeter Brune 
145bf7f4e0aSPeter Brune #undef __FUNCT__
146bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchDestroy"
147bf7f4e0aSPeter Brune PetscErrorCode LineSearchDestroy(LineSearch * linesearch) {
148bf7f4e0aSPeter Brune   PetscErrorCode ierr;
149bf7f4e0aSPeter Brune   PetscFunctionBegin;
150bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
151bf7f4e0aSPeter Brune   PetscValidHeaderSpecific((*linesearch),LineSearch_CLASSID,1);
152bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
153bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
154bf7f4e0aSPeter Brune   ierr = LineSearchReset(*linesearch);
155bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
156bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
157bf7f4e0aSPeter Brune   }
158bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
159e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
160bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
161bf7f4e0aSPeter Brune }
162bf7f4e0aSPeter Brune 
163bf7f4e0aSPeter Brune #undef __FUNCT__
164e7058c64SPeter Brune #define __FUNCT__ "LineSearchSetMonitor"
165bf7f4e0aSPeter Brune /*@C
166bf7f4e0aSPeter Brune    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search
167bf7f4e0aSPeter Brune 
168bf7f4e0aSPeter Brune    Input Parameters:
169bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
170bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
171bf7f4e0aSPeter Brune 
172bf7f4e0aSPeter Brune    Logically Collective on SNES
173bf7f4e0aSPeter Brune 
174bf7f4e0aSPeter Brune    Options Database:
175bf7f4e0aSPeter Brune .   -snes_ls_monitor
176bf7f4e0aSPeter Brune 
177bf7f4e0aSPeter Brune    Level: intermediate
178bf7f4e0aSPeter Brune 
179bf7f4e0aSPeter Brune 
180bf7f4e0aSPeter Brune .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
181bf7f4e0aSPeter Brune @*/
182bf7f4e0aSPeter Brune PetscErrorCode  LineSearchSetMonitor(LineSearch linesearch,PetscBool flg)
183bf7f4e0aSPeter Brune {
184bf7f4e0aSPeter Brune 
185bf7f4e0aSPeter Brune   PetscErrorCode ierr;
186bf7f4e0aSPeter Brune   PetscFunctionBegin;
187bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
188bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
189bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
190bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
191bf7f4e0aSPeter Brune   }
192bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
193bf7f4e0aSPeter Brune }
194bf7f4e0aSPeter Brune 
195bf7f4e0aSPeter Brune #undef __FUNCT__
196*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetMonitor"
197*6a388c36SPeter Brune 
198*6a388c36SPeter Brune PetscErrorCode  LineSearchGetMonitor(LineSearch linesearch, PetscViewer *monitor)
199*6a388c36SPeter Brune {
200*6a388c36SPeter Brune 
201*6a388c36SPeter Brune   PetscFunctionBegin;
202*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
203*6a388c36SPeter Brune   if (monitor) {
204*6a388c36SPeter Brune     PetscValidPointer(monitor, 2);
205*6a388c36SPeter Brune     *monitor = linesearch->monitor;
206*6a388c36SPeter Brune   }
207*6a388c36SPeter Brune   PetscFunctionReturn(0);
208*6a388c36SPeter Brune }
209*6a388c36SPeter Brune 
210*6a388c36SPeter Brune #undef __FUNCT__
211bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetFromOptions"
212bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetFromOptions(LineSearch linesearch) {
213bf7f4e0aSPeter Brune   PetscErrorCode ierr;
214bf7f4e0aSPeter Brune   const char     *deft = LINESEARCHBASIC;
215bf7f4e0aSPeter Brune   char           type[256];
216bf7f4e0aSPeter Brune   PetscBool      flg, set;
217bf7f4e0aSPeter Brune   PetscFunctionBegin;
218bf7f4e0aSPeter Brune   if (!LineSearchRegisterAllCalled) {ierr = LineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
219bf7f4e0aSPeter Brune 
220bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
221bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
222bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
223bf7f4e0aSPeter Brune   }
224bf7f4e0aSPeter Brune   ierr = PetscOptionsList("-linesearch_type","Line-search method","LineSearchSetType",LineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
225bf7f4e0aSPeter Brune   if (flg) {
226bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,type);CHKERRQ(ierr);
227bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
228bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,deft);CHKERRQ(ierr);
229bf7f4e0aSPeter Brune   }
230bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
231bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
232bf7f4e0aSPeter Brune   }
233bf7f4e0aSPeter Brune 
234bf7f4e0aSPeter Brune   ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESLineSearchSetMonitor",
235bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
236bf7f4e0aSPeter Brune   if (set) {ierr = LineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
237bf7f4e0aSPeter Brune 
238bf7f4e0aSPeter Brune   ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","LineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
239*6a388c36SPeter Brune   ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","LineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
240*6a388c36SPeter Brune   ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","LineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
241bf7f4e0aSPeter Brune   ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
242bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
243bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
244bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
245bf7f4e0aSPeter Brune }
246bf7f4e0aSPeter Brune 
247bf7f4e0aSPeter Brune #undef __FUNCT__
248bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchView"
249bf7f4e0aSPeter Brune PetscErrorCode LineSearchView(LineSearch linesearch) {
250bf7f4e0aSPeter Brune   PetscFunctionBegin;
251bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
252bf7f4e0aSPeter Brune }
253bf7f4e0aSPeter Brune 
254bf7f4e0aSPeter Brune #undef __FUNCT__
255bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetType"
256bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetType(LineSearch linesearch, const LineSearchType type)
257bf7f4e0aSPeter Brune {
258bf7f4e0aSPeter Brune 
259bf7f4e0aSPeter Brune   PetscErrorCode ierr,(*r)(LineSearch);
260bf7f4e0aSPeter Brune   PetscBool      match;
261bf7f4e0aSPeter Brune 
262bf7f4e0aSPeter Brune   PetscFunctionBegin;
263bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
264bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
265bf7f4e0aSPeter Brune 
266bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
267bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
268bf7f4e0aSPeter Brune 
269bf7f4e0aSPeter Brune   ierr =  PetscFListFind(LineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
270bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
271bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
272bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
273bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
274bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
275bf7f4e0aSPeter Brune   }
276bf7f4e0aSPeter Brune   /* Reinitialize function pointers in LineSearchOps structure */
277bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
278bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
279bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
280bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
281bf7f4e0aSPeter Brune 
282bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
283bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
284bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
285bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
286bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
287bf7f4e0aSPeter Brune   }
288bf7f4e0aSPeter Brune #endif
289bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
290bf7f4e0aSPeter Brune }
291bf7f4e0aSPeter Brune 
292bf7f4e0aSPeter Brune #undef __FUNCT__
293bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetSNES"
294bf7f4e0aSPeter Brune PetscErrorCode  LineSearchSetSNES(LineSearch linesearch, SNES snes){
295bf7f4e0aSPeter Brune   PetscFunctionBegin;
296bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
297bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
298bf7f4e0aSPeter Brune   linesearch->snes = snes;
299bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
300bf7f4e0aSPeter Brune }
301bf7f4e0aSPeter Brune 
302bf7f4e0aSPeter Brune #undef __FUNCT__
303bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetSNES"
304bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetSNES(LineSearch linesearch, SNES *snes){
305bf7f4e0aSPeter Brune   PetscFunctionBegin;
306*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
307*6a388c36SPeter Brune   PetscValidPointer(snes, 2);
308bf7f4e0aSPeter Brune   *snes = linesearch->snes;
309bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
310bf7f4e0aSPeter Brune }
311bf7f4e0aSPeter Brune 
312*6a388c36SPeter Brune #undef __FUNCT__
313*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetLambda"
314*6a388c36SPeter Brune PetscErrorCode  LineSearchGetLambda(LineSearch linesearch,PetscReal *lambda)
315*6a388c36SPeter Brune {
316*6a388c36SPeter Brune   PetscFunctionBegin;
317*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
318*6a388c36SPeter Brune   PetscValidPointer(lambda, 2);
319*6a388c36SPeter Brune   *lambda = linesearch->lambda;
320*6a388c36SPeter Brune   PetscFunctionReturn(0);
321*6a388c36SPeter Brune }
322*6a388c36SPeter Brune 
323*6a388c36SPeter Brune #undef __FUNCT__
324*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetLambda"
325*6a388c36SPeter Brune PetscErrorCode  LineSearchSetLambda(LineSearch linesearch, PetscReal lambda)
326*6a388c36SPeter Brune {
327*6a388c36SPeter Brune   PetscFunctionBegin;
328*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
329*6a388c36SPeter Brune   linesearch->lambda = lambda;
330*6a388c36SPeter Brune   PetscFunctionReturn(0);
331*6a388c36SPeter Brune }
332*6a388c36SPeter Brune 
333*6a388c36SPeter Brune #undef __FUNCT__
334*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetStepTolerance"
335*6a388c36SPeter Brune PetscErrorCode  LineSearchGetStepTolerance(LineSearch linesearch ,PetscReal *steptol)
336*6a388c36SPeter Brune {
337*6a388c36SPeter Brune   PetscFunctionBegin;
338*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
339*6a388c36SPeter Brune   PetscValidPointer(steptol, 2);
340*6a388c36SPeter Brune   *steptol = linesearch->steptol;
341*6a388c36SPeter Brune   PetscFunctionReturn(0);
342*6a388c36SPeter Brune }
343*6a388c36SPeter Brune 
344*6a388c36SPeter Brune #undef __FUNCT__
345*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetStepTolerance"
346*6a388c36SPeter Brune PetscErrorCode  LineSearchSetStepTolerance(LineSearch linesearch,PetscReal steptol)
347*6a388c36SPeter Brune {
348*6a388c36SPeter Brune   PetscFunctionBegin;
349*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
350*6a388c36SPeter Brune   linesearch->steptol = steptol;
351*6a388c36SPeter Brune   PetscFunctionReturn(0);
352*6a388c36SPeter Brune }
353*6a388c36SPeter Brune 
354*6a388c36SPeter Brune #undef __FUNCT__
355*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetDamping"
356*6a388c36SPeter Brune extern PetscErrorCode  LineSearchGetDamping(LineSearch linesearch,PetscReal *damping)
357*6a388c36SPeter Brune {
358*6a388c36SPeter Brune   PetscFunctionBegin;
359*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
360*6a388c36SPeter Brune   PetscValidPointer(damping, 2);
361*6a388c36SPeter Brune   *damping = linesearch->damping;
362*6a388c36SPeter Brune   PetscFunctionReturn(0);
363*6a388c36SPeter Brune }
364*6a388c36SPeter Brune 
365*6a388c36SPeter Brune #undef __FUNCT__
366*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetDamping"
367*6a388c36SPeter Brune extern PetscErrorCode  LineSearchSetDamping(LineSearch linesearch,PetscReal damping)
368*6a388c36SPeter Brune {
369*6a388c36SPeter Brune   PetscFunctionBegin;
370*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
371*6a388c36SPeter Brune   linesearch->damping = damping;
372*6a388c36SPeter Brune   PetscFunctionReturn(0);
373*6a388c36SPeter Brune }
374*6a388c36SPeter Brune 
375*6a388c36SPeter Brune #undef __FUNCT__
376*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetMaxStep"
377*6a388c36SPeter Brune extern PetscErrorCode  LineSearchGetMaxStep(LineSearch linesearch,PetscReal* maxstep)
378*6a388c36SPeter Brune {
379*6a388c36SPeter Brune   PetscFunctionBegin;
380*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
381*6a388c36SPeter Brune   PetscValidPointer(maxstep, 2);
382*6a388c36SPeter Brune   *maxstep = linesearch->maxstep;
383*6a388c36SPeter Brune   PetscFunctionReturn(0);
384*6a388c36SPeter Brune }
385*6a388c36SPeter Brune 
386*6a388c36SPeter Brune #undef __FUNCT__
387*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetMaxStep"
388*6a388c36SPeter Brune extern PetscErrorCode  LineSearchSetMaxStep(LineSearch linesearch, PetscReal maxstep)
389*6a388c36SPeter Brune {
390*6a388c36SPeter Brune   PetscFunctionBegin;
391*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
392*6a388c36SPeter Brune   linesearch->maxstep = maxstep;
393*6a388c36SPeter Brune   PetscFunctionReturn(0);
394*6a388c36SPeter Brune }
395*6a388c36SPeter Brune 
396*6a388c36SPeter Brune #undef __FUNCT__
397*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetMaxIts"
398*6a388c36SPeter Brune PetscErrorCode LineSearchGetMaxIts(LineSearch linesearch, PetscInt * max_its)
399*6a388c36SPeter Brune {
400*6a388c36SPeter Brune   PetscFunctionBegin;
401*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
402*6a388c36SPeter Brune   PetscValidPointer(max_its, 2);
403*6a388c36SPeter Brune   *max_its = linesearch->max_its;
404*6a388c36SPeter Brune   PetscFunctionReturn(0);
405*6a388c36SPeter Brune }
406*6a388c36SPeter Brune 
407*6a388c36SPeter Brune #undef __FUNCT__
408*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetMaxIts"
409*6a388c36SPeter Brune PetscErrorCode LineSearchSetMaxIts(LineSearch linesearch, PetscInt max_its)
410*6a388c36SPeter Brune {
411*6a388c36SPeter Brune   PetscFunctionBegin;
412*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
413*6a388c36SPeter Brune   linesearch->max_its = max_its;
414*6a388c36SPeter Brune   PetscFunctionReturn(0);
415*6a388c36SPeter Brune }
416bf7f4e0aSPeter Brune 
417bf7f4e0aSPeter Brune #undef __FUNCT__
418bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetNorms"
419bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetNorms(LineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
420bf7f4e0aSPeter Brune {
421bf7f4e0aSPeter Brune   PetscFunctionBegin;
422*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
423bf7f4e0aSPeter Brune   if (xnorm) {
424bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
425bf7f4e0aSPeter Brune   }
426bf7f4e0aSPeter Brune   if (fnorm) {
427bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
428bf7f4e0aSPeter Brune   }
429bf7f4e0aSPeter Brune   if (ynorm) {
430bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
431bf7f4e0aSPeter Brune   }
432bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
433bf7f4e0aSPeter Brune }
434bf7f4e0aSPeter Brune 
435e7058c64SPeter Brune #undef __FUNCT__
436*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetNorms"
437*6a388c36SPeter Brune PetscErrorCode  LineSearchSetNorms(LineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
438*6a388c36SPeter Brune {
439*6a388c36SPeter Brune   PetscFunctionBegin;
440*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
441*6a388c36SPeter Brune   if (xnorm) {
442*6a388c36SPeter Brune     linesearch->xnorm = xnorm;
443*6a388c36SPeter Brune   }
444*6a388c36SPeter Brune   if (fnorm) {
445*6a388c36SPeter Brune     linesearch->fnorm = fnorm;
446*6a388c36SPeter Brune   }
447*6a388c36SPeter Brune   if (ynorm) {
448*6a388c36SPeter Brune     linesearch->ynorm = ynorm;
449*6a388c36SPeter Brune   }
450*6a388c36SPeter Brune   PetscFunctionReturn(0);
451*6a388c36SPeter Brune }
452*6a388c36SPeter Brune 
453*6a388c36SPeter Brune #undef __FUNCT__
454*6a388c36SPeter Brune #define __FUNCT__ "LineSearchComputeNorms"
455*6a388c36SPeter Brune PetscErrorCode LineSearchComputeNorms(LineSearch linesearch)
456*6a388c36SPeter Brune {
457*6a388c36SPeter Brune   PetscErrorCode ierr;
458*6a388c36SPeter Brune   PetscFunctionBegin;
459*6a388c36SPeter Brune   if (linesearch->norms) {
460*6a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
461*6a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
462*6a388c36SPeter Brune     ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
463*6a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
464*6a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
465*6a388c36SPeter Brune     ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
466*6a388c36SPeter Brune   }
467*6a388c36SPeter Brune   PetscFunctionReturn(0);
468*6a388c36SPeter Brune }
469*6a388c36SPeter Brune 
470*6a388c36SPeter Brune #undef __FUNCT__
471*6a388c36SPeter Brune #define __FUNCT__ "LineSearchGetVecs"
472*6a388c36SPeter Brune PetscErrorCode LineSearchGetVecs(LineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
473*6a388c36SPeter Brune   PetscFunctionBegin;
474*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
475*6a388c36SPeter Brune   if (X) {
476*6a388c36SPeter Brune     PetscValidPointer(X, 2);
477*6a388c36SPeter Brune     *X = linesearch->vec_sol;
478*6a388c36SPeter Brune   }
479*6a388c36SPeter Brune   if (F) {
480*6a388c36SPeter Brune     PetscValidPointer(F, 3);
481*6a388c36SPeter Brune     *F = linesearch->vec_func;
482*6a388c36SPeter Brune   }
483*6a388c36SPeter Brune   if (Y) {
484*6a388c36SPeter Brune     PetscValidPointer(Y, 4);
485*6a388c36SPeter Brune     *Y = linesearch->vec_update;
486*6a388c36SPeter Brune   }
487*6a388c36SPeter Brune   if (W) {
488*6a388c36SPeter Brune     PetscValidPointer(W, 5);
489*6a388c36SPeter Brune     *W = linesearch->vec_sol_new;
490*6a388c36SPeter Brune   }
491*6a388c36SPeter Brune   if (G) {
492*6a388c36SPeter Brune     PetscValidPointer(G, 6);
493*6a388c36SPeter Brune     *G = linesearch->vec_func_new;
494*6a388c36SPeter Brune   }
495*6a388c36SPeter Brune 
496*6a388c36SPeter Brune   PetscFunctionReturn(0);
497*6a388c36SPeter Brune }
498*6a388c36SPeter Brune 
499*6a388c36SPeter Brune #undef __FUNCT__
500*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetVecs"
501*6a388c36SPeter Brune PetscErrorCode LineSearchSetVecs(LineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
502*6a388c36SPeter Brune   PetscFunctionBegin;
503*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
504*6a388c36SPeter Brune   if (X) {
505*6a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
506*6a388c36SPeter Brune     linesearch->vec_sol = X;
507*6a388c36SPeter Brune   }
508*6a388c36SPeter Brune   if (F) {
509*6a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
510*6a388c36SPeter Brune     linesearch->vec_func = F;
511*6a388c36SPeter Brune   }
512*6a388c36SPeter Brune   if (Y) {
513*6a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
514*6a388c36SPeter Brune     linesearch->vec_update = Y;
515*6a388c36SPeter Brune   }
516*6a388c36SPeter Brune   if (W) {
517*6a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
518*6a388c36SPeter Brune     linesearch->vec_sol_new = W;
519*6a388c36SPeter Brune   }
520*6a388c36SPeter Brune   if (G) {
521*6a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
522*6a388c36SPeter Brune     linesearch->vec_func_new = G;
523*6a388c36SPeter Brune   }
524*6a388c36SPeter Brune 
525*6a388c36SPeter Brune   PetscFunctionReturn(0);
526*6a388c36SPeter Brune }
527*6a388c36SPeter Brune 
528*6a388c36SPeter Brune #undef __FUNCT__
529e7058c64SPeter Brune #define __FUNCT__ "LineSearchAppendOptionsPrefix"
530e7058c64SPeter Brune /*@C
531e7058c64SPeter Brune    LineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
532e7058c64SPeter Brune    SNES options in the database.
533e7058c64SPeter Brune 
534e7058c64SPeter Brune    Logically Collective on SNES
535e7058c64SPeter Brune 
536e7058c64SPeter Brune    Input Parameters:
537e7058c64SPeter Brune +  snes - the SNES context
538e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
539e7058c64SPeter Brune 
540e7058c64SPeter Brune    Notes:
541e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
542e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
543e7058c64SPeter Brune 
544e7058c64SPeter Brune    Level: advanced
545e7058c64SPeter Brune 
546e7058c64SPeter Brune .keywords: SNES, append, options, prefix, database
547e7058c64SPeter Brune 
548e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
549e7058c64SPeter Brune @*/
550e7058c64SPeter Brune PetscErrorCode  LineSearchAppendOptionsPrefix(LineSearch linesearch,const char prefix[])
551e7058c64SPeter Brune {
552e7058c64SPeter Brune   PetscErrorCode ierr;
553e7058c64SPeter Brune 
554e7058c64SPeter Brune   PetscFunctionBegin;
555e7058c64SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
556e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
557e7058c64SPeter Brune   PetscFunctionReturn(0);
558e7058c64SPeter Brune }
559e7058c64SPeter Brune 
560e7058c64SPeter Brune #undef __FUNCT__
561e7058c64SPeter Brune #define __FUNCT__ "LineSearchGetOptionsPrefix"
562e7058c64SPeter Brune /*@C
563e7058c64SPeter Brune    LineSearchGetOptionsPrefix - Sets the prefix used for searching for all
564e7058c64SPeter Brune    SNES options in the database.
565e7058c64SPeter Brune 
566e7058c64SPeter Brune    Not Collective
567e7058c64SPeter Brune 
568e7058c64SPeter Brune    Input Parameter:
569e7058c64SPeter Brune .  snes - the SNES context
570e7058c64SPeter Brune 
571e7058c64SPeter Brune    Output Parameter:
572e7058c64SPeter Brune .  prefix - pointer to the prefix string used
573e7058c64SPeter Brune 
574e7058c64SPeter Brune    Notes: On the fortran side, the user should pass in a string 'prefix' of
575e7058c64SPeter Brune    sufficient length to hold the prefix.
576e7058c64SPeter Brune 
577e7058c64SPeter Brune    Level: advanced
578e7058c64SPeter Brune 
579e7058c64SPeter Brune .keywords: SNES, get, options, prefix, database
580e7058c64SPeter Brune 
581e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
582e7058c64SPeter Brune @*/
583e7058c64SPeter Brune PetscErrorCode  LineSearchGetOptionsPrefix(LineSearch linesearch,const char *prefix[])
584e7058c64SPeter Brune {
585e7058c64SPeter Brune   PetscErrorCode ierr;
586e7058c64SPeter Brune 
587e7058c64SPeter Brune   PetscFunctionBegin;
588e7058c64SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
589e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
590e7058c64SPeter Brune   PetscFunctionReturn(0);
591e7058c64SPeter Brune }
592bf7f4e0aSPeter Brune 
593bf7f4e0aSPeter Brune #undef __FUNCT__
594bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetWork"
595bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetWork(LineSearch linesearch, PetscInt nwork)
596bf7f4e0aSPeter Brune {
597bf7f4e0aSPeter Brune   PetscErrorCode ierr;
598bf7f4e0aSPeter Brune   PetscFunctionBegin;
599bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
600bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
601bf7f4e0aSPeter Brune   } else {
602bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
603bf7f4e0aSPeter Brune   }
604bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
605bf7f4e0aSPeter Brune }
606bf7f4e0aSPeter Brune 
607*6a388c36SPeter Brune 
608bf7f4e0aSPeter Brune #undef __FUNCT__
609bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetSuccess"
610bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetSuccess(LineSearch linesearch, PetscBool *success)
611bf7f4e0aSPeter Brune {
612bf7f4e0aSPeter Brune   PetscFunctionBegin;
613*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
614*6a388c36SPeter Brune   PetscValidPointer(success, 2);
615bf7f4e0aSPeter Brune   if (success) {
616bf7f4e0aSPeter Brune     *success = linesearch->success;
617bf7f4e0aSPeter Brune   }
618bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
619bf7f4e0aSPeter Brune }
620bf7f4e0aSPeter Brune 
621bf7f4e0aSPeter Brune #undef __FUNCT__
622*6a388c36SPeter Brune #define __FUNCT__ "LineSearchSetSuccess"
623*6a388c36SPeter Brune PetscErrorCode  LineSearchSetSuccess(LineSearch linesearch, PetscBool success)
624*6a388c36SPeter Brune {
625*6a388c36SPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
626*6a388c36SPeter Brune   PetscFunctionBegin;
627*6a388c36SPeter Brune   linesearch->success = success;
628*6a388c36SPeter Brune   PetscFunctionReturn(0);
629*6a388c36SPeter Brune }
630*6a388c36SPeter Brune 
631*6a388c36SPeter Brune #undef __FUNCT__
632bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchRegister"
633bf7f4e0aSPeter Brune /*@C
634bf7f4e0aSPeter Brune   LineSearchRegister - See LineSearchRegisterDynamic()
635bf7f4e0aSPeter Brune 
636bf7f4e0aSPeter Brune   Level: advanced
637bf7f4e0aSPeter Brune @*/
638bf7f4e0aSPeter Brune PetscErrorCode  LineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(LineSearch))
639bf7f4e0aSPeter Brune {
640bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
641bf7f4e0aSPeter Brune   PetscErrorCode ierr;
642bf7f4e0aSPeter Brune 
643bf7f4e0aSPeter Brune   PetscFunctionBegin;
644bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
645bf7f4e0aSPeter Brune   ierr = PetscFListAdd(&LineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
646bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
647bf7f4e0aSPeter Brune }
648