Actual source code: taosolver_bounds.c

  1: #include <petsc/private/taoimpl.h>

  3: /*@
  4:   TaoSetVariableBounds - Sets the upper and lower bounds

  6:   Logically collective on Tao

  8:   Input Parameters:
  9: + tao - the Tao context
 10: . XL  - vector of lower bounds
 11: - XU  - vector of upper bounds

 13:   Level: beginner

 15: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoGetVariableBounds()
 16: @*/
 17: PetscErrorCode TaoSetVariableBounds(Tao tao, Vec XL, Vec XU)
 18: {
 22:   PetscObjectReference((PetscObject)XL);
 23:   PetscObjectReference((PetscObject)XU);
 24:   VecDestroy(&tao->XL);
 25:   VecDestroy(&tao->XU);
 26:   tao->XL = XL;
 27:   tao->XU = XU;
 28:   tao->bounded = (PetscBool)(XL || XU);
 29:   return 0;
 30: }

 32: /*@C
 33:   TaoSetVariableBoundsRoutine - Sets a function to be used to compute variable bounds

 35:   Logically collective on Tao

 37:   Input Parameters:
 38: + tao - the Tao context
 39: . func - the bounds computation routine
 40: - ctx - [optional] user-defined context for private data for the bounds computation (may be NULL)

 42:   Calling sequence of func:
 43: $      func (Tao tao, Vec xl, Vec xu);

 45: + tao - the Tao
 46: . xl  - vector of lower bounds
 47: . xu  - vector of upper bounds
 48: - ctx - the (optional) user-defined function context

 50:   Level: beginner

 52: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetVariableBounds()

 54: Note: The func passed in to TaoSetVariableBoundsRoutine() takes
 55: precedence over any values set in TaoSetVariableBounds().

 57: @*/
 58: PetscErrorCode TaoSetVariableBoundsRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
 59: {
 61:   tao->user_boundsP = ctx;
 62:   tao->ops->computebounds = func;
 63:   tao->bounded = func ? PETSC_TRUE : PETSC_FALSE;
 64:   return 0;
 65: }

 67: /*@
 68:   TaoGetVariableBounds - Gets the upper and lower bounds vectors set with TaoSetVariableBounds

 70:   Not collective

 72:   Input Parameter:
 73: . tao - the Tao context

 75:   Output Parametrs:
 76: + XL  - vector of lower bounds
 77: - XU  - vector of upper bounds

 79:   Level: beginner

 81: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetVariableBounds()
 82: @*/
 83: PetscErrorCode TaoGetVariableBounds(Tao tao, Vec *XL, Vec *XU)
 84: {
 86:   if (XL) *XL = tao->XL;
 87:   if (XU) *XU = tao->XU;
 88:   return 0;
 89: }

 91: /*@C
 92:    TaoComputeVariableBounds - Compute the variable bounds using the
 93:    routine set by TaoSetVariableBoundsRoutine().

 95:    Collective on Tao

 97:    Input Parameter:
 98: .  tao - the Tao context

100:    Level: developer

102: .seealso: TaoSetVariableBoundsRoutine(), TaoSetVariableBounds()
103: @*/

105: PetscErrorCode TaoComputeVariableBounds(Tao tao)
106: {
108:   if (!tao->XL || !tao->XU) {
110:     VecDuplicate(tao->solution, &tao->XL);
111:     VecSet(tao->XL, PETSC_NINFINITY);
112:     VecDuplicate(tao->solution, &tao->XU);
113:     VecSet(tao->XU, PETSC_INFINITY);
114:   }
115:   if (tao->ops->computebounds) {
116:     PetscStackPush("Tao compute variable bounds");
117:     (*tao->ops->computebounds)(tao,tao->XL,tao->XU,tao->user_boundsP);
118:     PetscStackPop;
119:   }
120:   return 0;
121: }

123: /*@
124:   TaoSetInequalityBounds - Sets the upper and lower bounds

126:   Logically collective on Tao

128:   Input Parameters:
129: + tao - the Tao context
130: . IL  - vector of lower bounds
131: - IU  - vector of upper bounds

133:   Level: beginner

135: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoGetInequalityBounds()
136: @*/
137: PetscErrorCode TaoSetInequalityBounds(Tao tao, Vec IL, Vec IU)
138: {
142:   PetscObjectReference((PetscObject)IL);
143:   PetscObjectReference((PetscObject)IU);
144:   VecDestroy(&tao->IL);
145:   VecDestroy(&tao->IU);
146:   tao->IL = IL;
147:   tao->IU = IU;
148:   tao->ineq_doublesided = (PetscBool)(IL || IU);
149:   return 0;
150: }

152: /*@
153:   TaoGetInequalityBounds - Gets the upper and lower bounds set via TaoSetInequalityBounds

155:   Logically collective on Tao

157:   Input Parameter:
158: . tao - the Tao context

160:   Output Parameters:
161: + IL  - vector of lower bounds
162: - IU  - vector of upper bounds

164:   Level: beginner

166: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetInequalityBounds()
167: @*/
168: PetscErrorCode TaoGetInequalityBounds(Tao tao, Vec *IL, Vec *IU)
169: {
171:   if (IL) *IL = tao->IL;
172:   if (IU) *IU = tao->IU;
173:   return 0;
174: }

176: /*@C
177:    TaoComputeConstraints - Compute the variable bounds using the
178:    routine set by TaoSetConstraintsRoutine().

180:    Collective on Tao

182:    Input Parameters:
183: .  tao - the Tao context

185:    Level: developer

187: .seealso: TaoSetConstraintsRoutine(), TaoComputeJacobian()
188: @*/

190: PetscErrorCode TaoComputeConstraints(Tao tao, Vec X, Vec C)
191: {
198:   PetscLogEventBegin(TAO_ConstraintsEval,tao,X,C,NULL);
199:   PetscStackPush("Tao constraints evaluation routine");
200:   (*tao->ops->computeconstraints)(tao,X,C,tao->user_conP);
201:   PetscStackPop;
202:   PetscLogEventEnd(TAO_ConstraintsEval,tao,X,C,NULL);
203:   tao->nconstraints++;
204:   return 0;
205: }

207: /*@C
208:   TaoSetConstraintsRoutine - Sets a function to be used to compute constraints.  TAO only handles constraints under certain conditions, see manual for details

210:   Logically collective on Tao

212:   Input Parameters:
213: + tao - the Tao context
214: . c   - A vector that will be used to store constraint evaluation
215: . func - the bounds computation routine
216: - ctx - [optional] user-defined context for private data for the constraints computation (may be NULL)

218:   Calling sequence of func:
219: $      func (Tao tao, Vec x, Vec c, void *ctx);

221: + tao - the Tao
222: . x   - point to evaluate constraints
223: . c   - vector constraints evaluated at x
224: - ctx - the (optional) user-defined function context

226:   Level: intermediate

228: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetVariablevBounds()

230: @*/
231: PetscErrorCode TaoSetConstraintsRoutine(Tao tao, Vec c, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
232: {
235:   PetscObjectReference((PetscObject)c);
236:   VecDestroy(&tao->constraints);
237:   tao->constrained = func ? PETSC_TRUE : PETSC_FALSE;
238:   tao->constraints = c;
239:   tao->user_conP = ctx;
240:   tao->ops->computeconstraints = func;
241:   return 0;
242: }

244: /*@
245:   TaoComputeDualVariables - Computes the dual vectors corresponding to the bounds
246:   of the variables

248:   Collective on Tao

250:   Input Parameter:
251: . tao - the Tao context

253:   Output Parameters:
254: + DL - dual variable vector for the lower bounds
255: - DU - dual variable vector for the upper bounds

257:   Level: advanced

259:   Note:
260:   DL and DU should be created before calling this routine.  If calling
261:   this routine after using an unconstrained solver, DL and DU are set to all
262:   zeros.

264:   Level: advanced

266: .seealso: TaoComputeObjective(), TaoSetVariableBounds()
267: @*/
268: PetscErrorCode TaoComputeDualVariables(Tao tao, Vec DL, Vec DU)
269: {
275:   if (tao->ops->computedual) {
276:     (*tao->ops->computedual)(tao,DL,DU);
277:   } else {
278:     VecSet(DL,0.0);
279:     VecSet(DU,0.0);
280:   }
281:   return 0;
282: }

284: /*@
285:   TaoGetDualVariables - Gets pointers to the dual vectors

287:   Collective on Tao

289:   Input Parameter:
290: . tao - the Tao context

292:   Output Parameters:
293: + DE - dual variable vector for the lower bounds
294: - DI - dual variable vector for the upper bounds

296:   Level: advanced

298: .seealso: TaoComputeDualVariables()
299: @*/
300: PetscErrorCode TaoGetDualVariables(Tao tao, Vec *DE, Vec *DI)
301: {
303:   if (DE) *DE = tao->DE;
304:   if (DI) *DI = tao->DI;
305:   return 0;
306: }

308: /*@C
309:   TaoSetEqualityConstraintsRoutine - Sets a function to be used to compute constraints.  TAO only handles constraints under certain conditions, see manual for details

311:   Logically collective on Tao

313:   Input Parameters:
314: + tao - the Tao context
315: . ce   - A vector that will be used to store equality constraint evaluation
316: . func - the bounds computation routine
317: - ctx - [optional] user-defined context for private data for the equality constraints computation (may be NULL)

319:   Calling sequence of func:
320: $      func (Tao tao, Vec x, Vec ce, void *ctx);

322: + tao - the Tao
323: . x   - point to evaluate equality constraints
324: . ce   - vector of equality constraints evaluated at x
325: - ctx - the (optional) user-defined function context

327:   Level: intermediate

329: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetVariableBounds()

331: @*/
332: PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao tao, Vec ce, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
333: {
336:   PetscObjectReference((PetscObject)ce);
337:   VecDestroy(&tao->constraints_equality);
338:   tao->eq_constrained = func ? PETSC_TRUE : PETSC_FALSE;
339:   tao->constraints_equality = ce;
340:   tao->user_con_equalityP = ctx;
341:   tao->ops->computeequalityconstraints = func;
342:   return 0;
343: }

345: /*@C
346:   TaoSetInequalityConstraintsRoutine - Sets a function to be used to compute constraints.  TAO only handles constraints under certain conditions, see manual for details

348:   Logically collective on Tao

350:   Input Parameters:
351: + tao - the Tao context
352: . ci   - A vector that will be used to store inequality constraint evaluation
353: . func - the bounds computation routine
354: - ctx - [optional] user-defined context for private data for the inequality constraints computation (may be NULL)

356:   Calling sequence of func:
357: $      func (Tao tao, Vec x, Vec ci, void *ctx);

359: + tao - the Tao
360: . x   - point to evaluate inequality constraints
361: . ci   - vector of inequality constraints evaluated at x
362: - ctx - the (optional) user-defined function context

364:   Level: intermediate

366: .seealso: TaoSetObjective(), TaoSetHessian(), TaoSetObjectiveAndGradient(), TaoSetVariableBounds()

368: @*/
369: PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao tao, Vec ci, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
370: {
373:   PetscObjectReference((PetscObject)ci);
374:   VecDestroy(&tao->constraints_inequality);
375:   tao->constraints_inequality = ci;
376:   tao->ineq_constrained = func ? PETSC_TRUE : PETSC_FALSE;
377:   tao->user_con_inequalityP = ctx;
378:   tao->ops->computeinequalityconstraints = func;
379:   return 0;
380: }

382: /*@C
383:    TaoComputeEqualityConstraints - Compute the variable bounds using the
384:    routine set by TaoSetEqualityConstraintsRoutine().

386:    Collective on Tao

388:    Input Parameters:
389: .  tao - the Tao context

391:    Level: developer

393: .seealso: TaoSetEqualityConstraintsRoutine(), TaoComputeJacobianEquality()
394: @*/

396: PetscErrorCode TaoComputeEqualityConstraints(Tao tao, Vec X, Vec CE)
397: {
404:   PetscLogEventBegin(TAO_ConstraintsEval,tao,X,CE,NULL);
405:   PetscStackPush("Tao equality constraints evaluation routine");
406:   (*tao->ops->computeequalityconstraints)(tao,X,CE,tao->user_con_equalityP);
407:   PetscStackPop;
408:   PetscLogEventEnd(TAO_ConstraintsEval,tao,X,CE,NULL);
409:   tao->nconstraints++;
410:   return 0;
411: }

413: /*@C
414:    TaoComputeInequalityConstraints - Compute the variable bounds using the
415:    routine set by TaoSetInequalityConstraintsRoutine().

417:    Collective on Tao

419:    Input Parameters:
420: .  tao - the Tao context

422:    Level: developer

424: .seealso: TaoSetInequalityConstraintsRoutine(), TaoComputeJacobianInequality()
425: @*/

427: PetscErrorCode TaoComputeInequalityConstraints(Tao tao, Vec X, Vec CI)
428: {
435:   PetscLogEventBegin(TAO_ConstraintsEval,tao,X,CI,NULL);
436:   PetscStackPush("Tao inequality constraints evaluation routine");
437:   (*tao->ops->computeinequalityconstraints)(tao,X,CI,tao->user_con_inequalityP);
438:   PetscStackPop;
439:   PetscLogEventEnd(TAO_ConstraintsEval,tao,X,CI,NULL);
440:   tao->nconstraints++;
441:   return 0;
442: }