xref: /petsc/src/tao/interface/taosolver_fg.c (revision a7fac7c2b47dbcd8ece0cc2dfdfe0e63be1bb7b5)
1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "TaoSetInitialVector"
5 /*@
6   TaoSetInitialVector - Sets the initial guess for the solve
7 
8   Logically collective on Tao
9 
10   Input Parameters:
11 + tao - the Tao context
12 - x0  - the initial guess
13 
14   Level: beginner
15 .seealso: TaoCreate(), TaoSolve()
16 @*/
17 
18 PetscErrorCode TaoSetInitialVector(Tao tao, Vec x0)
19 {
20   PetscErrorCode ierr;
21 
22   PetscFunctionBegin;
23   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
24   if (x0) {
25     PetscValidHeaderSpecific(x0,VEC_CLASSID,2);
26     PetscObjectReference((PetscObject)x0);
27   }
28   ierr = VecDestroy(&tao->solution);CHKERRQ(ierr);
29   tao->solution = x0;
30   PetscFunctionReturn(0);
31 }
32 
33 #undef __FUNCT__
34 #define __FUNCT__ "TaoComputeGradient"
35 /*@
36   TaoComputeGradient - Computes the gradient of the objective function
37 
38   Collective on Tao
39 
40   Input Parameters:
41 + tao - the Tao context
42 - X - input vector
43 
44   Output Parameter:
45 . G - gradient vector
46 
47   Notes: TaoComputeGradient() is typically used within minimization implementations,
48   so most users would not generally call this routine themselves.
49 
50   Level: advanced
51 
52 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine()
53 @*/
54 PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
55 {
56   PetscErrorCode ierr;
57   PetscReal      dummy;
58 
59   PetscFunctionBegin;
60   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
61   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
62   PetscValidHeaderSpecific(G,VEC_CLASSID,2);
63   PetscCheckSameComm(tao,1,X,2);
64   PetscCheckSameComm(tao,1,G,3);
65   if (tao->ops->computegradient) {
66     ierr = PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
67     PetscStackPush("Tao user gradient evaluation routine");
68     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
69     PetscStackPop;
70     ierr = PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
71     tao->ngrads++;
72   } else if (tao->ops->computeobjectiveandgradient) {
73     ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr);
74     PetscStackPush("Tao user objective/gradient evaluation routine");
75     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP);CHKERRQ(ierr);
76     PetscStackPop;
77     ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr);
78     tao->nfuncgrads++;
79   }  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradientRoutine() has not been called");
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "TaoComputeObjective"
85 /*@
86   TaoComputeObjective - Computes the objective function value at a given point
87 
88   Collective on Tao
89 
90   Input Parameters:
91 + tao - the Tao context
92 - X - input vector
93 
94   Output Parameter:
95 . f - Objective value at X
96 
97   Notes: TaoComputeObjective() is typically used within minimization implementations,
98   so most users would not generally call this routine themselves.
99 
100   Level: advanced
101 
102 .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
103 @*/
104 PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
105 {
106   PetscErrorCode ierr;
107   Vec            temp;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
111   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
112   PetscCheckSameComm(tao,1,X,2);
113   if (tao->ops->computeobjective) {
114     ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
115     PetscStackPush("Tao user objective evaluation routine");
116     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
117     PetscStackPop;
118     ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
119     tao->nfuncs++;
120   } else if (tao->ops->computeobjectiveandgradient) {
121     ierr = PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n");CHKERRQ(ierr);
122     ierr = VecDuplicate(X,&temp);CHKERRQ(ierr);
123     ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,NULL,NULL);CHKERRQ(ierr);
124     PetscStackPush("Tao user objective/gradient evaluation routine");
125     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP);CHKERRQ(ierr);
126     PetscStackPop;
127     ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,NULL,NULL);CHKERRQ(ierr);
128     ierr = VecDestroy(&temp);CHKERRQ(ierr);
129     tao->nfuncgrads++;
130   }  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() has not been called");
131   ierr = PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "TaoComputeObjectiveAndGradient"
137 /*@
138   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point
139 
140   Collective on Tao
141 
142   Input Parameters:
143 + tao - the Tao context
144 - X - input vector
145 
146   Output Parameter:
147 + f - Objective value at X
148 - g - Gradient vector at X
149 
150   Notes: TaoComputeObjectiveAndGradient() is typically used within minimization implementations,
151   so most users would not generally call this routine themselves.
152 
153   Level: advanced
154 
155 .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
156 @*/
157 PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
158 {
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
163   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
164   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
165   PetscCheckSameComm(tao,1,X,2);
166   PetscCheckSameComm(tao,1,G,4);
167   if (tao->ops->computeobjectiveandgradient) {
168     ierr = PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr);
169     PetscStackPush("Tao user objective/gradient evaluation routine");
170     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP);CHKERRQ(ierr);
171     PetscStackPop;
172     if (tao->ops->computegradient == TaoDefaultComputeGradient) {
173       /* Overwrite gradient with finite difference gradient */
174       ierr = TaoDefaultComputeGradient(tao,X,G,tao->user_objgradP);CHKERRQ(ierr);
175     }
176     ierr = PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);CHKERRQ(ierr);
177     tao->nfuncgrads++;
178   } else if (tao->ops->computeobjective && tao->ops->computegradient) {
179     ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
180     PetscStackPush("Tao user objective evaluation routine");
181     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
182     PetscStackPop;
183     ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
184     tao->nfuncs++;
185     ierr = PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
186     PetscStackPush("Tao user gradient evaluation routine");
187     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
188     PetscStackPop;
189     ierr = PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
190     tao->ngrads++;
191   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() or TaoSetGradientRoutine() not set");
192   ierr = PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "TaoSetObjectiveRoutine"
198 /*@C
199   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization
200 
201   Logically collective on Tao
202 
203   Input Parameter:
204 + tao - the Tao context
205 . func - the objective function
206 - ctx - [optional] user-defined context for private data for the function evaluation
207         routine (may be NULL)
208 
209   Calling sequence of func:
210 $      func (Tao tao, Vec x, PetscReal *f, void *ctx);
211 
212 + x - input vector
213 . f - function value
214 - ctx - [optional] user-defined function context
215 
216   Level: beginner
217 
218 .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
219 @*/
220 PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
221 {
222   PetscFunctionBegin;
223   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
224   tao->user_objP = ctx;
225   tao->ops->computeobjective = func;
226   PetscFunctionReturn(0);
227 }
228 
229 #undef __FUNCT__
230 #define __FUNCT__ "TaoSetSeparableObjectiveRoutine"
231 /*@C
232   TaoSetSeparableObjectiveRoutine - Sets the function evaluation routine for least-square applications
233 
234   Logically collective on Tao
235 
236   Input Parameter:
237 + tao - the Tao context
238 . func - the objective function evaluation routine
239 - ctx - [optional] user-defined context for private data for the function evaluation
240         routine (may be NULL)
241 
242   Calling sequence of func:
243 $      func (Tao tao, Vec x, Vec f, void *ctx);
244 
245 + x - input vector
246 . f - function value vector
247 - ctx - [optional] user-defined function context
248 
249   Level: beginner
250 
251 .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
252 @*/
253 PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec sepobj, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
254 {
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
257   PetscValidHeaderSpecific(sepobj, VEC_CLASSID,2);
258   tao->user_sepobjP = ctx;
259   tao->sep_objective = sepobj;
260   tao->ops->computeseparableobjective = func;
261   PetscFunctionReturn(0);
262 }
263 
264 #undef __FUNCT__
265 #define __FUNCT__ "TaoComputeSeparableObjective"
266 /*@
267   TaoComputeSeparableObjective - Computes a separable objective function vector at a given point (for least-square applications)
268 
269   Collective on Tao
270 
271   Input Parameters:
272 + tao - the Tao context
273 - X - input vector
274 
275   Output Parameter:
276 . f - Objective vector at X
277 
278   Notes: TaoComputeSeparableObjective() is typically used within minimization implementations,
279   so most users would not generally call this routine themselves.
280 
281   Level: advanced
282 
283 .seealso: TaoSetSeparableObjectiveRoutine()
284 @*/
285 PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F)
286 {
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
291   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
292   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
293   PetscCheckSameComm(tao,1,X,2);
294   PetscCheckSameComm(tao,1,F,3);
295   if (tao->ops->computeseparableobjective) {
296     ierr = PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
297     PetscStackPush("Tao user separable objective evaluation routine");
298     ierr = (*tao->ops->computeseparableobjective)(tao,X,F,tao->user_sepobjP);CHKERRQ(ierr);
299     PetscStackPop;
300     ierr = PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
301     tao->nfuncs++;
302   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetSeparableObjectiveRoutine() has not been called");
303   ierr = PetscInfo(tao,"TAO separable function evaluation.\n");CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "TaoSetGradientRoutine"
309 /*@C
310   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization
311 
312   Logically collective on Tao
313 
314   Input Parameter:
315 + tao - the Tao context
316 . func - the gradient function
317 - ctx - [optional] user-defined context for private data for the gradient evaluation
318         routine (may be NULL)
319 
320   Calling sequence of func:
321 $      func (Tao tao, Vec x, Vec g, void *ctx);
322 
323 + x - input vector
324 . g - gradient value (output)
325 - ctx - [optional] user-defined function context
326 
327   Level: beginner
328 
329 .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
330 @*/
331 PetscErrorCode TaoSetGradientRoutine(Tao tao,  PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
332 {
333   PetscFunctionBegin;
334   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
335   tao->user_gradP = ctx;
336   tao->ops->computegradient = func;
337   PetscFunctionReturn(0);
338 }
339 
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "TaoSetObjectiveAndGradientRoutine"
343 /*@C
344   TaoSetObjectiveAndGradientRoutine - Sets a combined objective function and gradient evaluation routine for minimization
345 
346   Logically collective on Tao
347 
348   Input Parameter:
349 + tao - the Tao context
350 . func - the gradient function
351 - ctx - [optional] user-defined context for private data for the gradient evaluation
352         routine (may be NULL)
353 
354   Calling sequence of func:
355 $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
356 
357 + x - input vector
358 . f - objective value (output)
359 . g - gradient value (output)
360 - ctx - [optional] user-defined function context
361 
362   Level: beginner
363 
364 .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
365 @*/
366 PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx)
367 {
368   PetscFunctionBegin;
369   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
370   tao->user_objgradP = ctx;
371   tao->ops->computeobjectiveandgradient = func;
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "TaoIsObjectiveDefined"
377 /*@
378   TaoIsObjectiveDefined -- Checks to see if the user has
379   declared an objective-only routine.  Useful for determining when
380   it is appropriate to call TaoComputeObjective() or
381   TaoComputeObjectiveAndGradient()
382 
383   Collective on Tao
384 
385   Input Parameter:
386 + tao - the Tao context
387 - ctx - PETSC_TRUE if objective function routine is set by user,
388         PETSC_FALSE otherwise
389   Level: developer
390 
391 .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined()
392 @*/
393 PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
394 {
395   PetscFunctionBegin;
396   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
397   if (tao->ops->computeobjective == 0) *flg = PETSC_FALSE;
398   else *flg = PETSC_TRUE;
399   PetscFunctionReturn(0);
400 }
401 
402 #undef __FUNCT__
403 #define __FUNCT__ "TaoIsGradientDefined"
404 /*@
405   TaoIsGradientDefined -- Checks to see if the user has
406   declared an objective-only routine.  Useful for determining when
407   it is appropriate to call TaoComputeGradient() or
408   TaoComputeGradientAndGradient()
409 
410   Not Collective
411 
412   Input Parameter:
413 + tao - the Tao context
414 - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
415   Level: developer
416 
417 .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined()
418 @*/
419 PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
420 {
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
423   if (tao->ops->computegradient == 0) *flg = PETSC_FALSE;
424   else *flg = PETSC_TRUE;
425   PetscFunctionReturn(0);
426 }
427 
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "TaoIsObjectiveAndGradientDefined"
431 /*@
432   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
433   declared a joint objective/gradient routine.  Useful for determining when
434   it is appropriate to call TaoComputeObjective() or
435   TaoComputeObjectiveAndGradient()
436 
437   Not Collective
438 
439   Input Parameter:
440 + tao - the Tao context
441 - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise
442   Level: developer
443 
444 .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined()
445 @*/
446 PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
447 {
448   PetscFunctionBegin;
449   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
450   if (tao->ops->computeobjectiveandgradient == 0) *flg = PETSC_FALSE;
451   else *flg = PETSC_TRUE;
452   PetscFunctionReturn(0);
453 }
454 
455 
456 
457