Actual source code: taosolver_fg.c

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

 81: /*@
 82:   TaoComputeObjective - Computes the objective function value at a given point

 84:   Collective on Tao

 86:   Input Parameters:
 87: + tao - the Tao context
 88: - X - input vector

 90:   Output Parameter:
 91: . f - Objective value at X

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

 96:   Level: advanced

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

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

133: /*@
134:   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point

136:   Collective on Tao

138:   Input Parameters:
139: + tao - the Tao context
140: - X - input vector

142:   Output Parameter:
143: + f - Objective value at X
144: - g - Gradient vector at X

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

149:   Level: advanced

151: .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
152: @*/
153: PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
154: {

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

195: /*@C
196:   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization

198:   Logically collective on Tao

200:   Input Parameter:
201: + tao - the Tao context
202: . func - the objective function
203: - ctx - [optional] user-defined context for private data for the function evaluation
204:         routine (may be NULL)

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

209: + x - input vector
210: . f - function value
211: - ctx - [optional] user-defined function context

213:   Level: beginner

215: .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
216: @*/
217: PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
218: {
221:   tao->user_objP = ctx;
222:   tao->ops->computeobjective = func;
223:   return(0);
224: }

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

229:   Logically collective on Tao

231:   Input Parameter:
232: + tao - the Tao context
233: . func - the objective function evaluation routine
234: - ctx - [optional] user-defined context for private data for the function evaluation
235:         routine (may be NULL)

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

240: + x - input vector
241: . f - function value vector
242: - ctx - [optional] user-defined function context

244:   Level: beginner

246: .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
247: @*/
248: PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec sepobj, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
249: {
253:   tao->user_sepobjP = ctx;
254:   tao->sep_objective = sepobj;
255:   tao->ops->computeseparableobjective = func;
256:   return(0);
257: }

259: /*@
260:   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.

262:   Collective on Tao

264:   Input Parameters:
265: + tao - the Tao context
266: . sigma_v - vector of weights (diagonal terms only)
267: . n       - the number of weights (if using off-diagonal)
268: . rows    - index list of rows for sigma_w
269: . cols    - index list of columns for sigma_w
270: - vals - array of weights



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

276:   Level: intermediate

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

316:   Collective on Tao

318:   Input Parameters:
319: + tao - the Tao context
320: - X - input vector

322:   Output Parameter:
323: . f - Objective vector at X

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

328:   Level: advanced

330: .seealso: TaoSetSeparableObjectiveRoutine()
331: @*/
332: PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F)
333: {

342:   if (tao->ops->computeseparableobjective) {
343:     PetscLogEventBegin(Tao_ObjectiveEval,tao,X,NULL,NULL);
344:     PetscStackPush("Tao user separable objective evaluation routine");
345:     (*tao->ops->computeseparableobjective)(tao,X,F,tao->user_sepobjP);
346:     PetscStackPop;
347:     PetscLogEventEnd(Tao_ObjectiveEval,tao,X,NULL,NULL);
348:     tao->nfuncs++;
349:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetSeparableObjectiveRoutine() has not been called");
350:   PetscInfo(tao,"TAO separable function evaluation.\n");
351:   return(0);
352: }

354: /*@C
355:   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization

357:   Logically collective on Tao

359:   Input Parameter:
360: + tao - the Tao context
361: . func - the gradient function
362: - ctx - [optional] user-defined context for private data for the gradient evaluation
363:         routine (may be NULL)

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

368: + x - input vector
369: . g - gradient value (output)
370: - ctx - [optional] user-defined function context

372:   Level: beginner

374: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
375: @*/
376: PetscErrorCode TaoSetGradientRoutine(Tao tao,  PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
377: {
380:   tao->user_gradP = ctx;
381:   tao->ops->computegradient = func;
382:   return(0);
383: }

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

388:   Logically collective on Tao

390:   Input Parameter:
391: + tao - the Tao context
392: . func - the gradient function
393: - ctx - [optional] user-defined context for private data for the gradient evaluation
394:         routine (may be NULL)

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

399: + x - input vector
400: . f - objective value (output)
401: . g - gradient value (output)
402: - ctx - [optional] user-defined function context

404:   Level: beginner

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

417: /*@
418:   TaoIsObjectiveDefined -- Checks to see if the user has
419:   declared an objective-only routine.  Useful for determining when
420:   it is appropriate to call TaoComputeObjective() or
421:   TaoComputeObjectiveAndGradient()

423:   Collective on Tao

425:   Input Parameter:
426: + tao - the Tao context
427: - ctx - PETSC_TRUE if objective function routine is set by user,
428:         PETSC_FALSE otherwise
429:   Level: developer

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

442: /*@
443:   TaoIsGradientDefined -- Checks to see if the user has
444:   declared an objective-only routine.  Useful for determining when
445:   it is appropriate to call TaoComputeGradient() or
446:   TaoComputeGradientAndGradient()

448:   Not Collective

450:   Input Parameter:
451: + tao - the Tao context
452: - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
453:   Level: developer

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

466: /*@
467:   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
468:   declared a joint objective/gradient routine.  Useful for determining when
469:   it is appropriate to call TaoComputeObjective() or
470:   TaoComputeObjectiveAndGradient()

472:   Not Collective

474:   Input Parameter:
475: + tao - the Tao context
476: - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise
477:   Level: developer

479: .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined()
480: @*/
481: PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
482: {
485:   if (tao->ops->computeobjectiveandgradient == 0) *flg = PETSC_FALSE;
486:   else *flg = PETSC_TRUE;
487:   return(0);
488: }