Actual source code: taosolver_bounds.c

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

  5: /*@
  6:   TaoSetVariableBounds - Sets the upper and lower bounds

  8:   Logically collective on Tao

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

 15:   Level: beginner

 17: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
 18: @*/

 20: PetscErrorCode TaoSetVariableBounds(Tao tao, Vec XL, Vec XU)
 21: {

 26:   if (XL) {
 28:     PetscObjectReference((PetscObject)XL);
 29:   }
 30:   if (XU) {
 32:     PetscObjectReference((PetscObject)XU);
 33:   }
 34:   VecDestroy(&tao->XL);
 35:   VecDestroy(&tao->XU);
 36:   tao->XL = XL;
 37:   tao->XU = XU;
 38:   return(0);
 39: }

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

 46:   Logically collective on Tao

 48:   Input Parameters:
 49: + tao - the Tao context
 50: . func - the bounds computation routine
 51: - ctx - [optional] user-defined context for private data for the bounds computation (may be NULL)

 53:   Calling sequence of func:
 54: $      func (Tao tao, Vec xl, Vec xu);

 56: + tao - the Tao
 57: . xl  - vector of lower bounds
 58: . xu  - vector of upper bounds
 59: - ctx - the (optional) user-defined function context

 61:   Level: beginner

 63: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine(), TaoSetVariableBounds()

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

 68: @*/
 69: PetscErrorCode TaoSetVariableBoundsRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
 70: {
 73:   tao->user_boundsP = ctx;
 74:   tao->ops->computebounds = func;
 75:   return(0);
 76: }

 80: PetscErrorCode TaoGetVariableBounds(Tao tao, Vec *XL, Vec *XU)
 81: {
 84:   if (XL) {
 85:     *XL=tao->XL;
 86:   }
 87:   if (XU) {
 88:     *XU=tao->XU;
 89:   }
 90:   return(0);
 91: }

 95: /*@C
 96:    TaoComputeVariableBounds - Compute the variable bounds using the
 97:    routine set by TaoSetVariableBoundsRoutine().

 99:    Collective on Tao

101:    Input Parameters:
102: .  tao - the Tao context

104:    Level: developer

106: .seealso: TaoSetVariableBoundsRoutine(), TaoSetVariableBounds()
107: @*/

109: PetscErrorCode TaoComputeVariableBounds(Tao tao)
110: {

115:   if (!tao->ops->computebounds) return(0);
116:   if (!tao->XL || !tao->XU) {
117:     if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeVariableBounds");
118:     VecDuplicate(tao->solution, &tao->XL);
119:     VecSet(tao->XL, PETSC_NINFINITY);
120:     VecDuplicate(tao->solution, &tao->XU);
121:     VecSet(tao->XU, PETSC_INFINITY);
122:   }
123:   PetscStackPush("Tao compute variable bounds");
124:   (*tao->ops->computebounds)(tao,tao->XL,tao->XU,tao->user_boundsP);
125:   PetscStackPop;
126:   return(0);
127: }

131: /*@
132:   TaoSetInequalityBounds - Sets the upper and lower bounds

134:   Logically collective on Tao

136:   Input Parameters:
137: + tao - the Tao context
138: . IL  - vector of lower bounds
139: - IU  - vector of upper bounds

141:   Level: beginner

143: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
144: @*/

146: PetscErrorCode TaoSetInequalityBounds(Tao tao, Vec IL, Vec IU)
147: {

152:   if (IL) {
154:     PetscObjectReference((PetscObject)IL);
155:   }
156:   if (IU) {
158:     PetscObjectReference((PetscObject)IU);
159:   }
160:   VecDestroy(&tao->IL);
161:   VecDestroy(&tao->IU);
162:   tao->IL = IL;
163:   tao->IU = IU;
164:   return(0);
165: }


170: PetscErrorCode TaoGetInequalityBounds(Tao tao, Vec *IL, Vec *IU)
171: {
174:   if (IL) {
175:     *IL=tao->IL;
176:   }
177:   if (IU) {
178:     *IU=tao->IU;
179:   }
180:   return(0);
181: }

185: /*@C
186:    TaoComputeConstraints - Compute the variable bounds using the
187:    routine set by TaoSetConstraintsRoutine().

189:    Collective on Tao

191:    Input Parameters:
192: .  tao - the Tao context

194:    Level: developer

196: .seealso: TaoSetConstraintsRoutine(), TaoComputeJacobian()
197: @*/

199: PetscErrorCode TaoComputeConstraints(Tao tao, Vec X, Vec C)
200: {


210:   if (!tao->ops->computeconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetConstraintsRoutine() has not been called");
211:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeConstraints");
212:   PetscLogEventBegin(Tao_ConstraintsEval,tao,X,C,NULL);
213:   PetscStackPush("Tao constraints evaluation routine");
214:   (*tao->ops->computeconstraints)(tao,X,C,tao->user_conP);
215:   PetscStackPop;
216:   PetscLogEventEnd(Tao_ConstraintsEval,tao,X,C,NULL);
217:   tao->nconstraints++;
218:   return(0);
219: }

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

226:   Logically collective on Tao

228:   Input Parameters:
229: + tao - the Tao context
230: . c   - A vector that will be used to store constraint evaluation
231: . func - the bounds computation routine
232: - ctx - [optional] user-defined context for private data for the constraints computation (may be NULL)

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

237: + tao - the Tao
238: . x   - point to evaluate constraints
239: . c   - vector constraints evaluated at x
240: - ctx - the (optional) user-defined function context

242:   Level: intermediate

244: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine(), TaoSetVariablevBounds()

246: @*/
247: PetscErrorCode TaoSetConstraintsRoutine(Tao tao, Vec c, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
248: {
251:     tao->constraints = c;
252:     tao->user_conP = ctx;
253:     tao->ops->computeconstraints = func;
254:     return(0);
255: }

259: /*@
260:   TaoComputeDualVariables - Computes the dual vectors corresponding to the bounds
261:   of the variables

263:   Collective on Tao

265:   Input Parameters:
266: . tao - the Tao context

268:   Output Parameter:
269: + DL - dual variable vector for the lower bounds
270: - DU - dual variable vector for the upper bounds

272:   Level: advanced

274:   Note:
275:   DL and DU should be created before calling this routine.  If calling
276:   this routine after using an unconstrained solver, DL and DU are set to all
277:   zeros.

279:   Level: advanced

281: .seealso: TaoComputeObjective(), TaoSetVariableBounds()
282: @*/
283: PetscErrorCode TaoComputeDualVariables(Tao tao, Vec DL, Vec DU)
284: {
292:   if (tao->ops->computedual) {
293:     (*tao->ops->computedual)(tao,DL,DU);
294:   }  else {
295:     VecSet(DL,0.0);
296:     VecSet(DU,0.0);
297:   }
298:   return(0);
299: }

303: /*@
304:   TaoGetDualVariables - Gets pointers to the dual vectors

306:   Collective on Tao

308:   Input Parameters:
309: . tao - the Tao context

311:   Output Parameter:
312: + DE - dual variable vector for the lower bounds
313: - DI - dual variable vector for the upper bounds

315:   Level: advanced

317: .seealso: TaoComputeDualVariables()
318: @*/
319: PetscErrorCode TaoGetDualVariables(Tao tao, Vec *DE, Vec *DI)
320: {
323:   if (DE) {
324:     *DE = tao->DE;
325:   }
326:   if (DI) {
327:     *DI = tao->DI;
328:   }
329:   return(0);
330: }

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

337:   Logically collective on Tao

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

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

348: + tao - the Tao
349: . x   - point to evaluate equality constraints
350: . ce   - vector of equality constraints evaluated at x
351: - ctx - the (optional) user-defined function context

353:   Level: intermediate

355: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine(), TaoSetVariableBounds()

357: @*/
358: PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao tao, Vec ce, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
359: {

364:   if (ce) {
366:     PetscObjectReference((PetscObject)ce);
367:   }
368:   VecDestroy(&tao->constraints_equality);

370:   tao->constraints_equality = ce;
371:   tao->user_con_equalityP = ctx;
372:   tao->ops->computeequalityconstraints = func;
373:   return(0);
374: }


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

382:   Logically collective on Tao

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

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

393: + tao - the Tao
394: . x   - point to evaluate inequality constraints
395: . ci   - vector of inequality constraints evaluated at x
396: - ctx - the (optional) user-defined function context

398:   Level: intermediate

400: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine(), TaoSetVariableBounds()

402: @*/
403: PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao tao, Vec ci, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
404: {

409:   if (ci) {
411:     PetscObjectReference((PetscObject)ci);
412:   }
413:   VecDestroy(&tao->constraints_inequality);
414:   tao->constraints_inequality = ci;

416:   tao->user_con_inequalityP = ctx;
417:   tao->ops->computeinequalityconstraints = func;
418:   return(0);
419: }


424: /*@C
425:    TaoComputeEqualityConstraints - Compute the variable bounds using the
426:    routine set by TaoSetEqualityConstraintsRoutine().

428:    Collective on Tao

430:    Input Parameters:
431: .  tao - the Tao context

433:    Level: developer

435: .seealso: TaoSetEqualityConstraintsRoutine(), TaoComputeJacobianEquality()
436: @*/

438: PetscErrorCode TaoComputeEqualityConstraints(Tao tao, Vec X, Vec CE)
439: {


449:   if (!tao->ops->computeequalityconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetEqualityConstraintsRoutine() has not been called");
450:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeEqualityConstraints");
451:   PetscLogEventBegin(Tao_ConstraintsEval,tao,X,CE,NULL);
452:   PetscStackPush("Tao equality constraints evaluation routine");
453:   (*tao->ops->computeequalityconstraints)(tao,X,CE,tao->user_con_equalityP);
454:   PetscStackPop;
455:   PetscLogEventEnd(Tao_ConstraintsEval,tao,X,CE,NULL);
456:   tao->nconstraints++;
457:   return(0);
458: }


463: /*@C
464:    TaoComputeInequalityConstraints - Compute the variable bounds using the
465:    routine set by TaoSetInequalityConstraintsRoutine().

467:    Collective on Tao

469:    Input Parameters:
470: .  tao - the Tao context

472:    Level: developer

474: .seealso: TaoSetInequalityConstraintsRoutine(), TaoComputeJacobianInequality()
475: @*/

477: PetscErrorCode TaoComputeInequalityConstraints(Tao tao, Vec X, Vec CI)
478: {


488:   if (!tao->ops->computeinequalityconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInequalityConstraintsRoutine() has not been called");
489:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeInequalityConstraints");
490:   PetscLogEventBegin(Tao_ConstraintsEval,tao,X,CI,NULL);
491:   PetscStackPush("Tao inequality constraints evaluation routine");
492:   (*tao->ops->computeinequalityconstraints)(tao,X,CI,tao->user_con_inequalityP);
493:   PetscStackPop;
494:   PetscLogEventEnd(Tao_ConstraintsEval,tao,X,CI,NULL);
495:   tao->nconstraints++;
496:   return(0);
497: }