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