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