Actual source code: taosolver_fg.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/

  5: /*@
  6:   TaoSetInitialVector - Sets the initial guess for the solve

  8:   Logically collective on Tao

 10:   Input Parameters:
 11: + tao - the Tao context
 12: - x0  - the initial guess

 14:   Level: beginner
 15: .seealso: TaoCreate(), TaoSolve()
 16: @*/

 18: PetscErrorCode TaoSetInitialVector(Tao tao, Vec x0)
 19: {

 24:   if (x0) {
 26:     PetscObjectReference((PetscObject)x0);
 27:   }
 28:   VecDestroy(&tao->solution);
 29:   tao->solution = x0;
 30:   return(0);
 31: }

 35: /*@
 36:   TaoComputeGradient - Computes the gradient of the objective function

 38:   Collective on Tao

 40:   Input Parameters:
 41: + tao - the Tao context
 42: - X - input vector

 44:   Output Parameter:
 45: . G - gradient vector

 47:   Notes: TaoComputeGradient() is typically used within minimization implementations,
 48:   so most users would not generally call this routine themselves.

 50:   Level: advanced

 52: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine()
 53: @*/
 54: PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
 55: {
 57:   PetscReal      dummy;

 65:   if (tao->ops->computegradient) {
 66:     PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);
 67:     PetscStackPush("Tao user gradient evaluation routine");
 68:     (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);
 69:     PetscStackPop;
 70:     PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);
 71:     tao->ngrads++;
 72:   } else if (tao->ops->computeobjectiveandgradient) {
 73:     PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);
 74:     PetscStackPush("Tao user objective/gradient evaluation routine");
 75:     (*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP);
 76:     PetscStackPop;
 77:     PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);
 78:     tao->nfuncgrads++;
 79:   }  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradientRoutine() has not been called");
 80:   return(0);
 81: }

 85: /*@
 86:   TaoComputeObjective - Computes the objective function value at a given point

 88:   Collective on Tao

 90:   Input Parameters:
 91: + tao - the Tao context
 92: - X - input vector

 94:   Output Parameter:
 95: . f - Objective value at X

 97:   Notes: TaoComputeObjective() is typically used within minimization implementations,
 98:   so most users would not generally call this routine themselves.

100:   Level: advanced

102: .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
103: @*/
104: PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
105: {
107:   Vec            temp;

113:   if (tao->ops->computeobjective) {
114:     PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);
115:     PetscStackPush("Tao user objective evaluation routine");
116:     (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);
117:     PetscStackPop;
118:     PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);
119:     tao->nfuncs++;
120:   } else if (tao->ops->computeobjectiveandgradient) {
121:     PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n");
122:     VecDuplicate(X,&temp);
123:     PetscLogEventBegin(Tao_ObjGradientEval,tao,X,NULL,NULL);
124:     PetscStackPush("Tao user objective/gradient evaluation routine");
125:     (*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP);
126:     PetscStackPop;
127:     PetscLogEventEnd(Tao_ObjGradientEval,tao,X,NULL,NULL);
128:     VecDestroy(&temp);
129:     tao->nfuncgrads++;
130:   }  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() has not been called");
131:   PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));
132:   return(0);
133: }

137: /*@
138:   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point

140:   Collective on Tao

142:   Input Parameters:
143: + tao - the Tao context
144: - X - input vector

146:   Output Parameter:
147: + f - Objective value at X
148: - g - Gradient vector at X

150:   Notes: TaoComputeObjectiveAndGradient() is typically used within minimization implementations,
151:   so most users would not generally call this routine themselves.

153:   Level: advanced

155: .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
156: @*/
157: PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
158: {

167:   if (tao->ops->computeobjectiveandgradient) {
168:     PetscLogEventBegin(Tao_ObjGradientEval,tao,X,G,NULL);
169:     PetscStackPush("Tao user objective/gradient evaluation routine");
170:     (*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP);
171:     PetscStackPop;
172:     if (tao->ops->computegradient == TaoDefaultComputeGradient) {
173:       /* Overwrite gradient with finite difference gradient */
174:       TaoDefaultComputeGradient(tao,X,G,tao->user_objgradP);
175:     }
176:     PetscLogEventEnd(Tao_ObjGradientEval,tao,X,G,NULL);
177:     tao->nfuncgrads++;
178:   } else if (tao->ops->computeobjective && tao->ops->computegradient) {
179:     PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);
180:     PetscStackPush("Tao user objective evaluation routine");
181:     (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);
182:     PetscStackPop;
183:     PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);
184:     tao->nfuncs++;
185:     PetscLogEventBegin(Tao_GradientEval,tao,X,G,NULL);
186:     PetscStackPush("Tao user gradient evaluation routine");
187:     (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);
188:     PetscStackPop;
189:     PetscLogEventEnd(Tao_GradientEval,tao,X,G,NULL);
190:     tao->ngrads++;
191:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() or TaoSetGradientRoutine() not set");
192:   PetscInfo1(tao,"TAO Function evaluation: %14.12e\n",(double)(*f));
193:   return(0);
194: }

198: /*@C
199:   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization

201:   Logically collective on Tao

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)

209:   Calling sequence of func:
210: $      func (Tao tao, Vec x, PetscReal *f, void *ctx);

212: + x - input vector
213: . f - function value
214: - ctx - [optional] user-defined function context

216:   Level: beginner

218: .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
219: @*/
220: PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
221: {
224:   tao->user_objP = ctx;
225:   tao->ops->computeobjective = func;
226:   return(0);
227: }

231: /*@C
232:   TaoSetSeparableObjectiveRoutine - Sets the function evaluation routine for least-square applications

234:   Logically collective on Tao

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)

242:   Calling sequence of func:
243: $      func (Tao tao, Vec x, Vec f, void *ctx);

245: + x - input vector
246: . f - function value vector
247: - ctx - [optional] user-defined function context

249:   Level: beginner

251: .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
252: @*/
253: PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec sepobj, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
254: {
258:   tao->user_sepobjP = ctx;
259:   tao->sep_objective = sepobj;
260:   tao->ops->computeseparableobjective = func;
261:   return(0);
262: }

266: /*@
267:   TaoSetSeparableObjectiveWeights - Give weights for the separable objective values. A vector can be used if only diagonal terms are used, otherwise a matrix can be give. If this function is not used, or if sigma_v and sigma_w are both NULL, then the default identity matrix will be used for weights.

269:   Collective on Tao

271:   Input Parameters:
272: + tao - the Tao context
273: . sigma_v - vector of weights (diagonal terms only)
274: . n       - the number of weights (if using off-diagonal)
275: . rows    - index list of rows for sigma_w
276: . cols    - index list of columns for sigma_w
277: - vals - array of weights



281:   Note: Either sigma_v or sigma_w (or both) should be NULL

283:   Level: intermediate

285: .seealso: TaoSetSeparableObjectiveRoutine()
286: @*/
287: PetscErrorCode TaoSetSeparableObjectiveWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
288: {
290:   PetscInt       i;
293:   VecDestroy(&tao->sep_weights_v);
294:   tao->sep_weights_v=sigma_v;
295:   if (sigma_v) {
296:     PetscObjectReference((PetscObject)sigma_v);
297:   }
298:   if (vals) {
299:     if (tao->sep_weights_n) {
300:       PetscFree(tao->sep_weights_rows);
301:       PetscFree(tao->sep_weights_cols);
302:       PetscFree(tao->sep_weights_w);
303:     }
304:     PetscMalloc1(n,&tao->sep_weights_rows);
305:     PetscMalloc1(n,&tao->sep_weights_cols);
306:     PetscMalloc1(n,&tao->sep_weights_w);
307:     tao->sep_weights_n=n;
308:     for (i=0;i<n;i++) {
309:       tao->sep_weights_rows[i]=rows[i];
310:       tao->sep_weights_cols[i]=cols[i];
311:       tao->sep_weights_w[i]=vals[i];
312:     }
313:   } else {
314:     tao->sep_weights_n=0;
315:     tao->sep_weights_rows=0;
316:     tao->sep_weights_cols=0;
317:   }
318:   return(0);
319: }
322: /*@
323:   TaoComputeSeparableObjective - Computes a separable objective function vector at a given point (for least-square applications)

325:   Collective on Tao

327:   Input Parameters:
328: + tao - the Tao context
329: - X - input vector

331:   Output Parameter:
332: . f - Objective vector at X

334:   Notes: TaoComputeSeparableObjective() is typically used within minimization implementations,
335:   so most users would not generally call this routine themselves.

337:   Level: advanced

339: .seealso: TaoSetSeparableObjectiveRoutine()
340: @*/
341: PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F)
342: {

351:   if (tao->ops->computeseparableobjective) {
352:     PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);
353:     PetscStackPush("Tao user separable objective evaluation routine");
354:     (*tao->ops->computeseparableobjective)(tao,X,F,tao->user_sepobjP);
355:     PetscStackPop;
356:     PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);
357:     tao->nfuncs++;
358:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetSeparableObjectiveRoutine() has not been called");
359:   PetscInfo(tao,"TAO separable function evaluation.\n");
360:   return(0);
361: }

365: /*@C
366:   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization

368:   Logically collective on Tao

370:   Input Parameter:
371: + tao - the Tao context
372: . func - the gradient function
373: - ctx - [optional] user-defined context for private data for the gradient evaluation
374:         routine (may be NULL)

376:   Calling sequence of func:
377: $      func (Tao tao, Vec x, Vec g, void *ctx);

379: + x - input vector
380: . g - gradient value (output)
381: - ctx - [optional] user-defined function context

383:   Level: beginner

385: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
386: @*/
387: PetscErrorCode TaoSetGradientRoutine(Tao tao,  PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
388: {
391:   tao->user_gradP = ctx;
392:   tao->ops->computegradient = func;
393:   return(0);
394: }


399: /*@C
400:   TaoSetObjectiveAndGradientRoutine - Sets a combined objective function and gradient evaluation routine for minimization

402:   Logically collective on Tao

404:   Input Parameter:
405: + tao - the Tao context
406: . func - the gradient function
407: - ctx - [optional] user-defined context for private data for the gradient evaluation
408:         routine (may be NULL)

410:   Calling sequence of func:
411: $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);

413: + x - input vector
414: . f - objective value (output)
415: . g - gradient value (output)
416: - ctx - [optional] user-defined function context

418:   Level: beginner

420: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
421: @*/
422: PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx)
423: {
426:   tao->user_objgradP = ctx;
427:   tao->ops->computeobjectiveandgradient = func;
428:   return(0);
429: }

433: /*@
434:   TaoIsObjectiveDefined -- Checks to see if the user has
435:   declared an objective-only routine.  Useful for determining when
436:   it is appropriate to call TaoComputeObjective() or
437:   TaoComputeObjectiveAndGradient()

439:   Collective on Tao

441:   Input Parameter:
442: + tao - the Tao context
443: - ctx - PETSC_TRUE if objective function routine is set by user,
444:         PETSC_FALSE otherwise
445:   Level: developer

447: .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined()
448: @*/
449: PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
450: {
453:   if (tao->ops->computeobjective == 0) *flg = PETSC_FALSE;
454:   else *flg = PETSC_TRUE;
455:   return(0);
456: }

460: /*@
461:   TaoIsGradientDefined -- Checks to see if the user has
462:   declared an objective-only routine.  Useful for determining when
463:   it is appropriate to call TaoComputeGradient() or
464:   TaoComputeGradientAndGradient()

466:   Not Collective

468:   Input Parameter:
469: + tao - the Tao context
470: - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
471:   Level: developer

473: .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined()
474: @*/
475: PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
476: {
479:   if (tao->ops->computegradient == 0) *flg = PETSC_FALSE;
480:   else *flg = PETSC_TRUE;
481:   return(0);
482: }


487: /*@
488:   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
489:   declared a joint objective/gradient routine.  Useful for determining when
490:   it is appropriate to call TaoComputeObjective() or
491:   TaoComputeObjectiveAndGradient()

493:   Not Collective

495:   Input Parameter:
496: + tao - the Tao context
497: - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise
498:   Level: developer

500: .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined()
501: @*/
502: PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
503: {
506:   if (tao->ops->computeobjectiveandgradient == 0) *flg = PETSC_FALSE;
507:   else *flg = PETSC_TRUE;
508:   return(0);
509: }