xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision bf7f4e0aa335b99ab1a0e3a31a22f19305ab3c54)
1*bf7f4e0aSPeter Brune #include <private/linesearchimpl.h> /*I "petsclinesearch.h" I*/
2*bf7f4e0aSPeter Brune 
3*bf7f4e0aSPeter Brune PetscBool  LineSearchRegisterAllCalled = PETSC_FALSE;
4*bf7f4e0aSPeter Brune PetscFList LineSearchList              = PETSC_NULL;
5*bf7f4e0aSPeter Brune 
6*bf7f4e0aSPeter Brune PetscClassId   LineSearch_CLASSID;
7*bf7f4e0aSPeter Brune PetscLogEvent  LineSearch_Apply;
8*bf7f4e0aSPeter Brune 
9*bf7f4e0aSPeter Brune #undef __FUNCT__
10*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchCreate"
11*bf7f4e0aSPeter Brune PetscErrorCode LineSearchCreate(MPI_Comm comm, LineSearch * outlinesearch) {
12*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
13*bf7f4e0aSPeter Brune   LineSearch     linesearch;
14*bf7f4e0aSPeter Brune   PetscFunctionBegin;
15*bf7f4e0aSPeter Brune   ierr = PetscHeaderCreate(linesearch, _p_LineSearch,struct _LineSearchOps,LineSearch_CLASSID, 0,
16*bf7f4e0aSPeter Brune                            "LineSearch","Line-search method","LineSearch",comm,LineSearchDestroy,LineSearchView);CHKERRQ(ierr);
17*bf7f4e0aSPeter Brune 
18*bf7f4e0aSPeter Brune   linesearch->ops->precheckstep = PETSC_NULL;
19*bf7f4e0aSPeter Brune   linesearch->ops->postcheckstep = PETSC_NULL;
20*bf7f4e0aSPeter Brune 
21*bf7f4e0aSPeter Brune   linesearch->lambda        = 1.0;
22*bf7f4e0aSPeter Brune   linesearch->fnorm         = 1.0;
23*bf7f4e0aSPeter Brune   linesearch->ynorm         = 1.0;
24*bf7f4e0aSPeter Brune   linesearch->xnorm         = 1.0;
25*bf7f4e0aSPeter Brune   linesearch->success       = PETSC_TRUE;
26*bf7f4e0aSPeter Brune   linesearch->norms         = PETSC_TRUE;
27*bf7f4e0aSPeter Brune   linesearch->keeplambda    = PETSC_FALSE;
28*bf7f4e0aSPeter Brune   linesearch->damping       = 1.0;
29*bf7f4e0aSPeter Brune   linesearch->maxstep       = 1e8;
30*bf7f4e0aSPeter Brune   linesearch->steptol       = 1e-12;
31*bf7f4e0aSPeter Brune   linesearch->precheckctx   = PETSC_NULL;
32*bf7f4e0aSPeter Brune   linesearch->postcheckctx  = PETSC_NULL;
33*bf7f4e0aSPeter Brune   linesearch->max_its       = 1;
34*bf7f4e0aSPeter Brune   linesearch->setupcalled   = PETSC_FALSE;
35*bf7f4e0aSPeter Brune   *outlinesearch            = linesearch;
36*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
37*bf7f4e0aSPeter Brune }
38*bf7f4e0aSPeter Brune 
39*bf7f4e0aSPeter Brune #undef __FUNCT__
40*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetUp"
41*bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetUp(LineSearch linesearch) {
42*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
43*bf7f4e0aSPeter Brune   PetscFunctionBegin;
44*bf7f4e0aSPeter Brune 
45*bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
46*bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,LINESEARCHBASIC);CHKERRQ(ierr);
47*bf7f4e0aSPeter Brune   }
48*bf7f4e0aSPeter Brune 
49*bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
50*bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
51*bf7f4e0aSPeter Brune     ierr = VecDuplicate(linesearch->vec_func, &linesearch->vec_func_new);CHKERRQ(ierr);
52*bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
53*bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
54*bf7f4e0aSPeter Brune     }
55*bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
56*bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
57*bf7f4e0aSPeter Brune   }
58*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
59*bf7f4e0aSPeter Brune }
60*bf7f4e0aSPeter Brune 
61*bf7f4e0aSPeter Brune #undef __FUNCT__
62*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchReset"
63*bf7f4e0aSPeter Brune PetscErrorCode LineSearchReset(LineSearch linesearch) {
64*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
65*bf7f4e0aSPeter Brune   PetscFunctionBegin;
66*bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
67*bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
68*bf7f4e0aSPeter Brune   }
69*bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
70*bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
71*bf7f4e0aSPeter Brune 
72*bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
73*bf7f4e0aSPeter Brune   linesearch->nwork = 0;
74*bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
75*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
76*bf7f4e0aSPeter Brune }
77*bf7f4e0aSPeter Brune 
78*bf7f4e0aSPeter Brune #undef __FUNCT__
79*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchPreCheck"
80*bf7f4e0aSPeter Brune PetscErrorCode LineSearchPreCheck(LineSearch linesearch, PetscBool * changed)
81*bf7f4e0aSPeter Brune {
82*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
83*bf7f4e0aSPeter Brune   PetscFunctionBegin;
84*bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
85*bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
86*bf7f4e0aSPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed);CHKERRQ(ierr);
87*bf7f4e0aSPeter Brune   }
88*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
89*bf7f4e0aSPeter Brune }
90*bf7f4e0aSPeter Brune 
91*bf7f4e0aSPeter Brune #undef __FUNCT__
92*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchPostCheck"
93*bf7f4e0aSPeter Brune PetscErrorCode LineSearchPostCheck(LineSearch linesearch, PetscBool * changed_W, PetscBool * changed_Y)
94*bf7f4e0aSPeter Brune {
95*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
96*bf7f4e0aSPeter Brune   PetscFunctionBegin;
97*bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
98*bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
99*bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
100*bf7f4e0aSPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_sol_new, linesearch->vec_update, changed_W, changed_Y);CHKERRQ(ierr);
101*bf7f4e0aSPeter Brune   }
102*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
103*bf7f4e0aSPeter Brune }
104*bf7f4e0aSPeter Brune 
105*bf7f4e0aSPeter Brune #undef __FUNCT__
106*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchApply"
107*bf7f4e0aSPeter Brune PetscErrorCode LineSearchApply(LineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
108*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
109*bf7f4e0aSPeter Brune   PetscFunctionBegin;
110*bf7f4e0aSPeter Brune 
111*bf7f4e0aSPeter Brune   /* check the pointers */
112*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
113*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
114*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
115*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
116*bf7f4e0aSPeter Brune 
117*bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
118*bf7f4e0aSPeter Brune 
119*bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
120*bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
121*bf7f4e0aSPeter Brune   linesearch->vec_func = F;
122*bf7f4e0aSPeter Brune 
123*bf7f4e0aSPeter Brune   ierr = LineSearchSetUp(linesearch);CHKERRQ(ierr);
124*bf7f4e0aSPeter Brune 
125*bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
126*bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
127*bf7f4e0aSPeter Brune 
128*bf7f4e0aSPeter Brune   if (fnorm) {
129*bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
130*bf7f4e0aSPeter Brune   } else {
131*bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
132*bf7f4e0aSPeter Brune   }
133*bf7f4e0aSPeter Brune 
134*bf7f4e0aSPeter Brune   ierr = PetscLogEventBegin(LineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
135*bf7f4e0aSPeter Brune 
136*bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
137*bf7f4e0aSPeter Brune 
138*bf7f4e0aSPeter Brune   ierr = PetscLogEventEnd(LineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
139*bf7f4e0aSPeter Brune 
140*bf7f4e0aSPeter Brune   if (fnorm)
141*bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
142*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
143*bf7f4e0aSPeter Brune }
144*bf7f4e0aSPeter Brune 
145*bf7f4e0aSPeter Brune #undef __FUNCT__
146*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchDestroy"
147*bf7f4e0aSPeter Brune PetscErrorCode LineSearchDestroy(LineSearch * linesearch) {
148*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
149*bf7f4e0aSPeter Brune   PetscFunctionBegin;
150*bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
151*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific((*linesearch),LineSearch_CLASSID,1);
152*bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
153*bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
154*bf7f4e0aSPeter Brune   ierr = LineSearchReset(*linesearch);
155*bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
156*bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
157*bf7f4e0aSPeter Brune   }
158*bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
159*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
160*bf7f4e0aSPeter Brune }
161*bf7f4e0aSPeter Brune 
162*bf7f4e0aSPeter Brune #undef __FUNCT__
163*bf7f4e0aSPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
164*bf7f4e0aSPeter Brune /*@C
165*bf7f4e0aSPeter Brune    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search
166*bf7f4e0aSPeter Brune 
167*bf7f4e0aSPeter Brune    Input Parameters:
168*bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
169*bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
170*bf7f4e0aSPeter Brune 
171*bf7f4e0aSPeter Brune    Logically Collective on SNES
172*bf7f4e0aSPeter Brune 
173*bf7f4e0aSPeter Brune    Options Database:
174*bf7f4e0aSPeter Brune .   -snes_ls_monitor
175*bf7f4e0aSPeter Brune 
176*bf7f4e0aSPeter Brune    Level: intermediate
177*bf7f4e0aSPeter Brune 
178*bf7f4e0aSPeter Brune 
179*bf7f4e0aSPeter Brune .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
180*bf7f4e0aSPeter Brune @*/
181*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchSetMonitor(LineSearch linesearch,PetscBool flg)
182*bf7f4e0aSPeter Brune {
183*bf7f4e0aSPeter Brune 
184*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
185*bf7f4e0aSPeter Brune   PetscFunctionBegin;
186*bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
187*bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
188*bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
189*bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
190*bf7f4e0aSPeter Brune   }
191*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
192*bf7f4e0aSPeter Brune }
193*bf7f4e0aSPeter Brune 
194*bf7f4e0aSPeter Brune #undef __FUNCT__
195*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetFromOptions"
196*bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetFromOptions(LineSearch linesearch) {
197*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
198*bf7f4e0aSPeter Brune   const char     *deft = LINESEARCHBASIC;
199*bf7f4e0aSPeter Brune   char           type[256];
200*bf7f4e0aSPeter Brune   PetscBool      flg, set;
201*bf7f4e0aSPeter Brune   PetscFunctionBegin;
202*bf7f4e0aSPeter Brune   if (!LineSearchRegisterAllCalled) {ierr = LineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
203*bf7f4e0aSPeter Brune 
204*bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
205*bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
206*bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
207*bf7f4e0aSPeter Brune   }
208*bf7f4e0aSPeter Brune   ierr = PetscOptionsList("-linesearch_type","Line-search method","LineSearchSetType",LineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
209*bf7f4e0aSPeter Brune   if (flg) {
210*bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,type);CHKERRQ(ierr);
211*bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
212*bf7f4e0aSPeter Brune     ierr = LineSearchSetType(linesearch,deft);CHKERRQ(ierr);
213*bf7f4e0aSPeter Brune   }
214*bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
215*bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
216*bf7f4e0aSPeter Brune   }
217*bf7f4e0aSPeter Brune 
218*bf7f4e0aSPeter Brune     ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESLineSearchSetMonitor",
219*bf7f4e0aSPeter Brune                             linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
220*bf7f4e0aSPeter Brune     if (set) {ierr = LineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
221*bf7f4e0aSPeter Brune 
222*bf7f4e0aSPeter Brune   ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","LineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
223*bf7f4e0aSPeter Brune   ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","LineSearchSetDamping",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
224*bf7f4e0aSPeter Brune   ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","LineSearchSetDamping",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
225*bf7f4e0aSPeter Brune   ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
226*bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
227*bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
228*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
229*bf7f4e0aSPeter Brune }
230*bf7f4e0aSPeter Brune 
231*bf7f4e0aSPeter Brune #undef __FUNCT__
232*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchView"
233*bf7f4e0aSPeter Brune PetscErrorCode LineSearchView(LineSearch linesearch) {
234*bf7f4e0aSPeter Brune   PetscFunctionBegin;
235*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
236*bf7f4e0aSPeter Brune }
237*bf7f4e0aSPeter Brune 
238*bf7f4e0aSPeter Brune #undef __FUNCT__
239*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetType"
240*bf7f4e0aSPeter Brune PetscErrorCode LineSearchSetType(LineSearch linesearch, const LineSearchType type)
241*bf7f4e0aSPeter Brune {
242*bf7f4e0aSPeter Brune 
243*bf7f4e0aSPeter Brune   PetscErrorCode ierr,(*r)(LineSearch);
244*bf7f4e0aSPeter Brune   PetscBool      match;
245*bf7f4e0aSPeter Brune 
246*bf7f4e0aSPeter Brune   PetscFunctionBegin;
247*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
248*bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
249*bf7f4e0aSPeter Brune 
250*bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
251*bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
252*bf7f4e0aSPeter Brune 
253*bf7f4e0aSPeter Brune   ierr =  PetscFListFind(LineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
254*bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
255*bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
256*bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
257*bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
258*bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
259*bf7f4e0aSPeter Brune   }
260*bf7f4e0aSPeter Brune   /* Reinitialize function pointers in LineSearchOps structure */
261*bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
262*bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
263*bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
264*bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
265*bf7f4e0aSPeter Brune 
266*bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
267*bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
268*bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
269*bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
270*bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
271*bf7f4e0aSPeter Brune   }
272*bf7f4e0aSPeter Brune #endif
273*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
274*bf7f4e0aSPeter Brune }
275*bf7f4e0aSPeter Brune 
276*bf7f4e0aSPeter Brune #undef __FUNCT__
277*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchSetSNES"
278*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchSetSNES(LineSearch linesearch, SNES snes){
279*bf7f4e0aSPeter Brune   PetscFunctionBegin;
280*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(linesearch,LineSearch_CLASSID,1);
281*bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
282*bf7f4e0aSPeter Brune   linesearch->snes = snes;
283*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
284*bf7f4e0aSPeter Brune }
285*bf7f4e0aSPeter Brune 
286*bf7f4e0aSPeter Brune #undef __FUNCT__
287*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetSNES"
288*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetSNES(LineSearch linesearch, SNES *snes){
289*bf7f4e0aSPeter Brune   PetscFunctionBegin;
290*bf7f4e0aSPeter Brune   *snes = linesearch->snes;
291*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
292*bf7f4e0aSPeter Brune }
293*bf7f4e0aSPeter Brune 
294*bf7f4e0aSPeter Brune 
295*bf7f4e0aSPeter Brune #undef __FUNCT__
296*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetNorms"
297*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetNorms(LineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
298*bf7f4e0aSPeter Brune {
299*bf7f4e0aSPeter Brune   PetscFunctionBegin;
300*bf7f4e0aSPeter Brune   if (xnorm) {
301*bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
302*bf7f4e0aSPeter Brune   }
303*bf7f4e0aSPeter Brune   if (fnorm) {
304*bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
305*bf7f4e0aSPeter Brune   }
306*bf7f4e0aSPeter Brune   if (ynorm) {
307*bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
308*bf7f4e0aSPeter Brune   }
309*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
310*bf7f4e0aSPeter Brune }
311*bf7f4e0aSPeter Brune 
312*bf7f4e0aSPeter Brune 
313*bf7f4e0aSPeter Brune #undef __FUNCT__
314*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetWork"
315*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetWork(LineSearch linesearch, PetscInt nwork)
316*bf7f4e0aSPeter Brune {
317*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
318*bf7f4e0aSPeter Brune   PetscFunctionBegin;
319*bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
320*bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
321*bf7f4e0aSPeter Brune   } else {
322*bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
323*bf7f4e0aSPeter Brune   }
324*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
325*bf7f4e0aSPeter Brune }
326*bf7f4e0aSPeter Brune 
327*bf7f4e0aSPeter Brune #undef __FUNCT__
328*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchGetSuccess"
329*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchGetSuccess(LineSearch linesearch, PetscBool *success)
330*bf7f4e0aSPeter Brune {
331*bf7f4e0aSPeter Brune   PetscFunctionBegin;
332*bf7f4e0aSPeter Brune   if (success) {
333*bf7f4e0aSPeter Brune     *success = linesearch->success;
334*bf7f4e0aSPeter Brune   }
335*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
336*bf7f4e0aSPeter Brune }
337*bf7f4e0aSPeter Brune 
338*bf7f4e0aSPeter Brune #undef __FUNCT__
339*bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchRegister"
340*bf7f4e0aSPeter Brune /*@C
341*bf7f4e0aSPeter Brune   LineSearchRegister - See LineSearchRegisterDynamic()
342*bf7f4e0aSPeter Brune 
343*bf7f4e0aSPeter Brune   Level: advanced
344*bf7f4e0aSPeter Brune @*/
345*bf7f4e0aSPeter Brune PetscErrorCode  LineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(LineSearch))
346*bf7f4e0aSPeter Brune {
347*bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
348*bf7f4e0aSPeter Brune   PetscErrorCode ierr;
349*bf7f4e0aSPeter Brune 
350*bf7f4e0aSPeter Brune   PetscFunctionBegin;
351*bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
352*bf7f4e0aSPeter Brune   ierr = PetscFListAdd(&LineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
353*bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
354*bf7f4e0aSPeter Brune }
355