Actual source code: taosolver_bounds.c

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

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

  6:   Logically Collective

  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: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetVariableBounds()`
 16: @*/
 17: PetscErrorCode TaoSetVariableBounds(Tao tao, Vec XL, Vec XU)
 18: {
 19:   PetscFunctionBegin;
 23:   PetscCall(PetscObjectReference((PetscObject)XL));
 24:   PetscCall(PetscObjectReference((PetscObject)XU));
 25:   PetscCall(VecDestroy(&tao->XL));
 26:   PetscCall(VecDestroy(&tao->XU));
 27:   tao->XL      = XL;
 28:   tao->XU      = XU;
 29:   tao->bounded = (PetscBool)(XL || XU);
 30:   PetscFunctionReturn(PETSC_SUCCESS);
 31: }

 33: /*@C
 34:   TaoSetVariableBoundsRoutine - Sets a function to be used to compute lower and upper variable bounds for the optimization

 36:   Logically Collective

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

 43:   Calling sequence of `func`:
 44: + tao - the `Tao` solver
 45: . xl  - vector of lower bounds
 46: . xu  - vector of upper bounds
 47: - ctx - the (optional) user-defined function context

 49:   Level: beginner

 51:   Note:
 52:   The func passed to `TaoSetVariableBoundsRoutine()` takes precedence over any values set in `TaoSetVariableBounds()`.

 54: .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetVariableBounds()`
 55: @*/
 56: PetscErrorCode TaoSetVariableBoundsRoutine(Tao tao, PetscErrorCode (*func)(Tao tao, Vec xl, Vec xu, void *ctx), void *ctx)
 57: {
 58:   PetscFunctionBegin;
 60:   tao->user_boundsP       = ctx;
 61:   tao->ops->computebounds = func;
 62:   tao->bounded            = func ? PETSC_TRUE : PETSC_FALSE;
 63:   PetscFunctionReturn(PETSC_SUCCESS);
 64: }

 66: /*@
 67:   TaoGetVariableBounds - Gets the upper and lower bounds vectors set with `TaoSetVariableBounds()`

 69:   Not Collective

 71:   Input Parameter:
 72: . tao - the `Tao` context

 74:   Output Parameters:
 75: + XL - vector of lower bounds
 76: - XU - vector of upper bounds

 78:   Level: beginner

 80: .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetVariableBounds()`
 81: @*/
 82: PetscErrorCode TaoGetVariableBounds(Tao tao, Vec *XL, Vec *XU)
 83: {
 84:   PetscFunctionBegin;
 86:   if (XL) *XL = tao->XL;
 87:   if (XU) *XU = tao->XU;
 88:   PetscFunctionReturn(PETSC_SUCCESS);
 89: }

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

 95:   Collective

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

100:   Level: developer

102: .seealso: [](ch_tao), `Tao`, `TaoSetVariableBoundsRoutine()`, `TaoSetVariableBounds()`
103: @*/
104: PetscErrorCode TaoComputeVariableBounds(Tao tao)
105: {
106:   PetscFunctionBegin;
108:   if (tao->ops->computebounds) {
109:     if (!tao->XL) {
110:       PetscCall(VecDuplicate(tao->solution, &tao->XL));
111:       PetscCall(VecSet(tao->XL, PETSC_NINFINITY));
112:     }
113:     if (!tao->XU) {
114:       PetscCall(VecDuplicate(tao->solution, &tao->XU));
115:       PetscCall(VecSet(tao->XU, PETSC_INFINITY));
116:     }
117:     PetscCallBack("Tao callback variable bounds", (*tao->ops->computebounds)(tao, tao->XL, tao->XU, tao->user_boundsP));
118:   }
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

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

125:   Logically Collective

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

132:   Level: beginner

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

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

155:   Logically Collective

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: [](ch_tao), `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetInequalityBounds()`
167: @*/
168: PetscErrorCode TaoGetInequalityBounds(Tao tao, Vec *IL, Vec *IU)
169: {
170:   PetscFunctionBegin;
172:   if (IL) *IL = tao->IL;
173:   if (IU) *IU = tao->IU;
174:   PetscFunctionReturn(PETSC_SUCCESS);
175: }

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

181:   Collective

183:   Input Parameters:
184: + tao - the `Tao` context
185: - X   - location to evaluate the constraints

187:   Output Parameter:
188: . C - the constraints

190:   Level: developer

192: .seealso: [](ch_tao), `Tao`, `TaoSetConstraintsRoutine()`, `TaoComputeJacobian()`
193: @*/
194: PetscErrorCode TaoComputeConstraints(Tao tao, Vec X, Vec C)
195: {
196:   PetscFunctionBegin;
200:   PetscCheckSameComm(tao, 1, X, 2);
201:   PetscCheckSameComm(tao, 1, C, 3);
202:   PetscCall(PetscLogEventBegin(TAO_ConstraintsEval, tao, X, C, NULL));
203:   PetscCallBack("Tao callback constraints", (*tao->ops->computeconstraints)(tao, X, C, tao->user_conP));
204:   PetscCall(PetscLogEventEnd(TAO_ConstraintsEval, tao, X, C, NULL));
205:   tao->nconstraints++;
206:   PetscFunctionReturn(PETSC_SUCCESS);
207: }

209: /*@C
210:   TaoSetConstraintsRoutine - Sets a function to be used to compute constraints.  Tao only handles constraints under certain conditions, see [](ch_tao) for details

212:   Logically Collective

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

220:   Calling sequence of `func`:
221: + tao - the `Tao` solver
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: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetVariablevBounds()`
229: @*/
230: PetscErrorCode TaoSetConstraintsRoutine(Tao tao, Vec c, PetscErrorCode (*func)(Tao tao, Vec x, Vec c, void *ctx), void *ctx)
231: {
232:   PetscFunctionBegin;
235:   PetscCall(PetscObjectReference((PetscObject)c));
236:   PetscCall(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:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

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

248:   Collective

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: .seealso: [](ch_tao), `Tao`, `TaoComputeObjective()`, `TaoSetVariableBounds()`
265: @*/
266: PetscErrorCode TaoComputeDualVariables(Tao tao, Vec DL, Vec DU)
267: {
268:   PetscFunctionBegin;
272:   PetscCheckSameComm(tao, 1, DL, 2);
273:   PetscCheckSameComm(tao, 1, DU, 3);
274:   if (tao->ops->computedual) {
275:     PetscUseTypeMethod(tao, computedual, DL, DU);
276:   } else {
277:     PetscCall(VecSet(DL, 0.0));
278:     PetscCall(VecSet(DU, 0.0));
279:   }
280:   PetscFunctionReturn(PETSC_SUCCESS);
281: }

283: /*@
284:   TaoGetDualVariables - Gets the dual vectors

286:   Collective

288:   Input Parameter:
289: . tao - the `Tao` context

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

295:   Level: advanced

297: .seealso: [](ch_tao), `Tao`, `TaoComputeDualVariables()`
298: @*/
299: PetscErrorCode TaoGetDualVariables(Tao tao, Vec *DE, Vec *DI)
300: {
301:   PetscFunctionBegin;
303:   if (DE) *DE = tao->DE;
304:   if (DI) *DI = tao->DI;
305:   PetscFunctionReturn(PETSC_SUCCESS);
306: }

308: /*@C
309:   TaoSetEqualityConstraintsRoutine - Sets a function to be used to compute constraints.  Tao only handles constraints under certain conditions, see [](ch_tao) for details

311:   Logically Collective

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: + tao - the `Tao` solver
321: . x   - point to evaluate equality constraints
322: . ce  - vector of equality constraints evaluated at x
323: - ctx - the (optional) user-defined function context

325:   Level: intermediate

327: .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetVariableBounds()`
328: @*/
329: PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao tao, Vec ce, PetscErrorCode (*func)(Tao tao, Vec x, Vec ce, void *ctx), void *ctx)
330: {
331:   PetscFunctionBegin;
334:   PetscCall(PetscObjectReference((PetscObject)ce));
335:   PetscCall(VecDestroy(&tao->constraints_equality));
336:   tao->eq_constrained                  = func ? PETSC_TRUE : PETSC_FALSE;
337:   tao->constraints_equality            = ce;
338:   tao->user_con_equalityP              = ctx;
339:   tao->ops->computeequalityconstraints = func;
340:   PetscFunctionReturn(PETSC_SUCCESS);
341: }

343: /*@C
344:   TaoSetInequalityConstraintsRoutine - Sets a function to be used to compute constraints.  Tao only handles constraints under certain conditions, see [](ch_tao) for details

346:   Logically Collective

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

354:   Calling sequence of `func`:
355: + tao - the `Tao` solver
356: . x   - point to evaluate inequality constraints
357: . ci  - vector of inequality constraints evaluated at x
358: - ctx - the (optional) user-defined function context

360:   Level: intermediate

362: .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetVariableBounds()`
363: @*/
364: PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao tao, Vec ci, PetscErrorCode (*func)(Tao tao, Vec x, Vec ci, void *ctx), void *ctx)
365: {
366:   PetscFunctionBegin;
369:   PetscCall(PetscObjectReference((PetscObject)ci));
370:   PetscCall(VecDestroy(&tao->constraints_inequality));
371:   tao->constraints_inequality            = ci;
372:   tao->ineq_constrained                  = func ? PETSC_TRUE : PETSC_FALSE;
373:   tao->user_con_inequalityP              = ctx;
374:   tao->ops->computeinequalityconstraints = func;
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: /*@C
379:   TaoComputeEqualityConstraints - Compute the variable bounds using the
380:   routine set by `TaoSetEqualityConstraintsRoutine()`.

382:   Collective

384:   Input Parameter:
385: . tao - the `Tao` context

387:   Output Parameters:
388: + X  - point the equality constraints were evaluated on
389: - CE - vector of equality constraints evaluated at X

391:   Level: developer

393: .seealso: [](ch_tao), `Tao`, `TaoSetEqualityConstraintsRoutine()`, `TaoComputeJacobianEquality()`, `TaoComputeInequalityConstraints()`
394: @*/
395: PetscErrorCode TaoComputeEqualityConstraints(Tao tao, Vec X, Vec CE)
396: {
397:   PetscFunctionBegin;
401:   PetscCheckSameComm(tao, 1, X, 2);
402:   PetscCheckSameComm(tao, 1, CE, 3);
403:   PetscCall(PetscLogEventBegin(TAO_ConstraintsEval, tao, X, CE, NULL));
404:   PetscCallBack("Tao callback equality constraints", (*tao->ops->computeequalityconstraints)(tao, X, CE, tao->user_con_equalityP));
405:   PetscCall(PetscLogEventEnd(TAO_ConstraintsEval, tao, X, CE, NULL));
406:   tao->nconstraints++;
407:   PetscFunctionReturn(PETSC_SUCCESS);
408: }

410: /*@C
411:   TaoComputeInequalityConstraints - Compute the variable bounds using the
412:   routine set by `TaoSetInequalityConstraintsRoutine()`.

414:   Collective

416:   Input Parameter:
417: . tao - the `Tao` context

419:   Output Parameters:
420: + X  - point the inequality constraints were evaluated on
421: - CI - vector of inequality constraints evaluated at X

423:   Level: developer

425: .seealso: [](ch_tao), `Tao`, `TaoSetInequalityConstraintsRoutine()`, `TaoComputeJacobianInequality()`, `TaoComputeEqualityConstraints()`
426: @*/
427: PetscErrorCode TaoComputeInequalityConstraints(Tao tao, Vec X, Vec CI)
428: {
429:   PetscFunctionBegin;
433:   PetscCheckSameComm(tao, 1, X, 2);
434:   PetscCheckSameComm(tao, 1, CI, 3);
435:   PetscCall(PetscLogEventBegin(TAO_ConstraintsEval, tao, X, CI, NULL));
436:   PetscCallBack("Tao callback inequality constraints", (*tao->ops->computeinequalityconstraints)(tao, X, CI, tao->user_con_inequalityP));
437:   PetscCall(PetscLogEventEnd(TAO_ConstraintsEval, tao, X, CI, NULL));
438:   tao->nconstraints++;
439:   PetscFunctionReturn(PETSC_SUCCESS);
440: }