Actual source code: taosolver_fg.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1:  #include <petsc/private/taoimpl.h>

  3: /*@
  4:   TaoSetInitialVector - Sets the initial guess for the solve

  6:   Logically collective on Tao

  8:   Input Parameters:
  9: + tao - the Tao context
 10: - x0  - the initial guess

 12:   Level: beginner
 13: .seealso: TaoCreate(), TaoSolve()
 14: @*/

 16: PetscErrorCode TaoSetInitialVector(Tao tao, Vec x0)
 17: {

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

 31: /*@
 32:   TaoComputeGradient - Computes the gradient of the objective function

 34:   Collective on Tao

 36:   Input Parameters:
 37: + tao - the Tao context
 38: - X - input vector

 40:   Output Parameter:
 41: . G - gradient vector

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

 46:   Level: advanced

 48: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine()
 49: @*/
 50: PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
 51: {
 53:   PetscReal      dummy;

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

 79: /*@
 80:   TaoComputeObjective - Computes the objective function value at a given point

 82:   Collective on Tao

 84:   Input Parameters:
 85: + tao - the Tao context
 86: - X - input vector

 88:   Output Parameter:
 89: . f - Objective value at X

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

 94:   Level: advanced

 96: .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
 97: @*/
 98: PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
 99: {
101:   Vec            temp;

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

129: /*@
130:   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point

132:   Collective on Tao

134:   Input Parameters:
135: + tao - the Tao context
136: - X - input vector

138:   Output Parameter:
139: + f - Objective value at X
140: - g - Gradient vector at X

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

145:   Level: advanced

147: .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
148: @*/
149: PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
150: {

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

188: /*@C
189:   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization

191:   Logically collective on Tao

193:   Input Parameter:
194: + tao - the Tao context
195: . func - the objective function
196: - ctx - [optional] user-defined context for private data for the function evaluation
197:         routine (may be NULL)

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

202: + x - input vector
203: . f - function value
204: - ctx - [optional] user-defined function context

206:   Level: beginner

208: .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
209: @*/
210: PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
211: {
214:   tao->user_objP = ctx;
215:   tao->ops->computeobjective = func;
216:   return(0);
217: }

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

222:   Logically collective on Tao

224:   Input Parameter:
225: + tao - the Tao context
226: . func - the objective function evaluation routine
227: - ctx - [optional] user-defined context for private data for the function evaluation
228:         routine (may be NULL)

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

233: + x - input vector
234: . f - function value vector
235: - ctx - [optional] user-defined function context

237:   Level: beginner

239: .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
240: @*/
241: PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec sepobj, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
242: {
246:   tao->user_sepobjP = ctx;
247:   tao->sep_objective = sepobj;
248:   tao->ops->computeseparableobjective = func;
249:   return(0);
250: }

252: /*@
253:   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.

255:   Collective on Tao

257:   Input Parameters:
258: + tao - the Tao context
259: . sigma_v - vector of weights (diagonal terms only)
260: . n       - the number of weights (if using off-diagonal)
261: . rows    - index list of rows for sigma_w
262: . cols    - index list of columns for sigma_w
263: - vals - array of weights



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

269:   Level: intermediate

271: .seealso: TaoSetSeparableObjectiveRoutine()
272: @*/
273: PetscErrorCode TaoSetSeparableObjectiveWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
274: {
276:   PetscInt       i;
279:   VecDestroy(&tao->sep_weights_v);
280:   tao->sep_weights_v=sigma_v;
281:   if (sigma_v) {
282:     PetscObjectReference((PetscObject)sigma_v);
283:   }
284:   if (vals) {
285:     if (tao->sep_weights_n) {
286:       PetscFree(tao->sep_weights_rows);
287:       PetscFree(tao->sep_weights_cols);
288:       PetscFree(tao->sep_weights_w);
289:     }
290:     PetscMalloc1(n,&tao->sep_weights_rows);
291:     PetscMalloc1(n,&tao->sep_weights_cols);
292:     PetscMalloc1(n,&tao->sep_weights_w);
293:     tao->sep_weights_n=n;
294:     for (i=0;i<n;i++) {
295:       tao->sep_weights_rows[i]=rows[i];
296:       tao->sep_weights_cols[i]=cols[i];
297:       tao->sep_weights_w[i]=vals[i];
298:     }
299:   } else {
300:     tao->sep_weights_n=0;
301:     tao->sep_weights_rows=0;
302:     tao->sep_weights_cols=0;
303:   }
304:   return(0);
305: }
306: /*@
307:   TaoComputeSeparableObjective - Computes a separable objective function vector at a given point (for least-square applications)

309:   Collective on Tao

311:   Input Parameters:
312: + tao - the Tao context
313: - X - input vector

315:   Output Parameter:
316: . f - Objective vector at X

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

321:   Level: advanced

323: .seealso: TaoSetSeparableObjectiveRoutine()
324: @*/
325: PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F)
326: {

335:   if (tao->ops->computeseparableobjective) {
336:     PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);
337:     PetscStackPush("Tao user separable objective evaluation routine");
338:     (*tao->ops->computeseparableobjective)(tao,X,F,tao->user_sepobjP);
339:     PetscStackPop;
340:     PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);
341:     tao->nfuncs++;
342:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetSeparableObjectiveRoutine() has not been called");
343:   PetscInfo(tao,"TAO separable function evaluation.\n");
344:   return(0);
345: }

347: /*@C
348:   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization

350:   Logically collective on Tao

352:   Input Parameter:
353: + tao - the Tao context
354: . func - the gradient function
355: - ctx - [optional] user-defined context for private data for the gradient evaluation
356:         routine (may be NULL)

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

361: + x - input vector
362: . g - gradient value (output)
363: - ctx - [optional] user-defined function context

365:   Level: beginner

367: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
368: @*/
369: PetscErrorCode TaoSetGradientRoutine(Tao tao,  PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
370: {
373:   tao->user_gradP = ctx;
374:   tao->ops->computegradient = func;
375:   return(0);
376: }


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

382:   Logically collective on Tao

384:   Input Parameter:
385: + tao - the Tao context
386: . func - the gradient function
387: - ctx - [optional] user-defined context for private data for the gradient evaluation
388:         routine (may be NULL)

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

393: + x - input vector
394: . f - objective value (output)
395: . g - gradient value (output)
396: - ctx - [optional] user-defined function context

398:   Level: beginner

400: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
401: @*/
402: PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx)
403: {
406:   tao->user_objgradP = ctx;
407:   tao->ops->computeobjectiveandgradient = func;
408:   return(0);
409: }

411: /*@
412:   TaoIsObjectiveDefined -- Checks to see if the user has
413:   declared an objective-only routine.  Useful for determining when
414:   it is appropriate to call TaoComputeObjective() or
415:   TaoComputeObjectiveAndGradient()

417:   Collective on Tao

419:   Input Parameter:
420: + tao - the Tao context
421: - ctx - PETSC_TRUE if objective function routine is set by user,
422:         PETSC_FALSE otherwise
423:   Level: developer

425: .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined()
426: @*/
427: PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
428: {
431:   if (tao->ops->computeobjective == 0) *flg = PETSC_FALSE;
432:   else *flg = PETSC_TRUE;
433:   return(0);
434: }

436: /*@
437:   TaoIsGradientDefined -- Checks to see if the user has
438:   declared an objective-only routine.  Useful for determining when
439:   it is appropriate to call TaoComputeGradient() or
440:   TaoComputeGradientAndGradient()

442:   Not Collective

444:   Input Parameter:
445: + tao - the Tao context
446: - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
447:   Level: developer

449: .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined()
450: @*/
451: PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
452: {
455:   if (tao->ops->computegradient == 0) *flg = PETSC_FALSE;
456:   else *flg = PETSC_TRUE;
457:   return(0);
458: }


461: /*@
462:   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
463:   declared a joint objective/gradient routine.  Useful for determining when
464:   it is appropriate to call TaoComputeObjective() or
465:   TaoComputeObjectiveAndGradient()

467:   Not Collective

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

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