Actual source code: taosolver_bounds.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  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: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
 16: @*/

 18: PetscErrorCode TaoSetVariableBounds(Tao tao, Vec XL, Vec XU)
 19: {

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

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

 42:   Logically collective on Tao

 44:   Input Parameters:
 45: + tao - the Tao context
 46: . func - the bounds computation routine
 47: - ctx - [optional] user-defined context for private data for the bounds computation (may be NULL)

 49:   Calling sequence of func:
 50: $      func (Tao tao, Vec xl, Vec xu);

 52: + tao - the Tao
 53: . xl  - vector of lower bounds
 54: . xu  - vector of upper bounds
 55: - ctx - the (optional) user-defined function context

 57:   Level: beginner

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

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

 64: @*/
 65: PetscErrorCode TaoSetVariableBoundsRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
 66: {
 69:   tao->user_boundsP = ctx;
 70:   tao->ops->computebounds = func;
 71:   return(0);
 72: }

 74: PetscErrorCode TaoGetVariableBounds(Tao tao, Vec *XL, Vec *XU)
 75: {
 78:   if (XL) {
 79:     *XL=tao->XL;
 80:   }
 81:   if (XU) {
 82:     *XU=tao->XU;
 83:   }
 84:   return(0);
 85: }

 87: /*@C
 88:    TaoComputeVariableBounds - Compute the variable bounds using the
 89:    routine set by TaoSetVariableBoundsRoutine().

 91:    Collective on Tao

 93:    Input Parameters:
 94: .  tao - the Tao context

 96:    Level: developer

 98: .seealso: TaoSetVariableBoundsRoutine(), TaoSetVariableBounds()
 99: @*/

101: PetscErrorCode TaoComputeVariableBounds(Tao tao)
102: {

107:   if (!tao->ops->computebounds) return(0);
108:   if (!tao->XL || !tao->XU) {
109:     if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeVariableBounds");
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:   PetscStackPush("Tao compute variable bounds");
116:   (*tao->ops->computebounds)(tao,tao->XL,tao->XU,tao->user_boundsP);
117:   PetscStackPop;
118:   return(0);
119: }

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

124:   Logically collective on Tao

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

131:   Level: beginner

133: .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
134: @*/

136: PetscErrorCode TaoSetInequalityBounds(Tao tao, Vec IL, Vec IU)
137: {

142:   if (IL) {
144:     PetscObjectReference((PetscObject)IL);
145:   }
146:   if (IU) {
148:     PetscObjectReference((PetscObject)IU);
149:   }
150:   VecDestroy(&tao->IL);
151:   VecDestroy(&tao->IU);
152:   tao->IL = IL;
153:   tao->IU = IU;
154:   return(0);
155: }


158: PetscErrorCode TaoGetInequalityBounds(Tao tao, Vec *IL, Vec *IU)
159: {
162:   if (IL) {
163:     *IL=tao->IL;
164:   }
165:   if (IU) {
166:     *IU=tao->IU;
167:   }
168:   return(0);
169: }

171: /*@C
172:    TaoComputeConstraints - Compute the variable bounds using the
173:    routine set by TaoSetConstraintsRoutine().

175:    Collective on Tao

177:    Input Parameters:
178: .  tao - the Tao context

180:    Level: developer

182: .seealso: TaoSetConstraintsRoutine(), TaoComputeJacobian()
183: @*/

185: PetscErrorCode TaoComputeConstraints(Tao tao, Vec X, Vec C)
186: {


196:   if (!tao->ops->computeconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetConstraintsRoutine() has not been called");
197:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeConstraints");
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: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine(), TaoSetVariablevBounds()

230: @*/
231: PetscErrorCode TaoSetConstraintsRoutine(Tao tao, Vec c, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
232: {
235:     tao->constraints = c;
236:     tao->user_conP = ctx;
237:     tao->ops->computeconstraints = func;
238:     return(0);
239: }

241: /*@
242:   TaoComputeDualVariables - Computes the dual vectors corresponding to the bounds
243:   of the variables

245:   Collective on Tao

247:   Input Parameters:
248: . tao - the Tao context

250:   Output Parameter:
251: + DL - dual variable vector for the lower bounds
252: - DU - dual variable vector for the upper bounds

254:   Level: advanced

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

261:   Level: advanced

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

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

286:   Collective on Tao

288:   Input Parameters:
289: . tao - the Tao context

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

295:   Level: advanced

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

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

315:   Logically collective on Tao

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

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

326: + tao - the Tao
327: . x   - point to evaluate equality constraints
328: . ce   - vector of equality constraints evaluated at x
329: - ctx - the (optional) user-defined function context

331:   Level: intermediate

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

335: @*/
336: PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao tao, Vec ce, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
337: {

342:   if (ce) {
344:     PetscObjectReference((PetscObject)ce);
345:   }
346:   VecDestroy(&tao->constraints_equality);

348:   tao->constraints_equality = ce;
349:   tao->user_con_equalityP = ctx;
350:   tao->ops->computeequalityconstraints = func;
351:   return(0);
352: }


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

358:   Logically collective on Tao

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

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

369: + tao - the Tao
370: . x   - point to evaluate inequality constraints
371: . ci   - vector of inequality constraints evaluated at x
372: - ctx - the (optional) user-defined function context

374:   Level: intermediate

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

378: @*/
379: PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao tao, Vec ci, PetscErrorCode (*func)(Tao, Vec, Vec, void*), void *ctx)
380: {

385:   if (ci) {
387:     PetscObjectReference((PetscObject)ci);
388:   }
389:   VecDestroy(&tao->constraints_inequality);
390:   tao->constraints_inequality = ci;

392:   tao->user_con_inequalityP = ctx;
393:   tao->ops->computeinequalityconstraints = func;
394:   return(0);
395: }


398: /*@C
399:    TaoComputeEqualityConstraints - Compute the variable bounds using the
400:    routine set by TaoSetEqualityConstraintsRoutine().

402:    Collective on Tao

404:    Input Parameters:
405: .  tao - the Tao context

407:    Level: developer

409: .seealso: TaoSetEqualityConstraintsRoutine(), TaoComputeJacobianEquality()
410: @*/

412: PetscErrorCode TaoComputeEqualityConstraints(Tao tao, Vec X, Vec CE)
413: {


423:   if (!tao->ops->computeequalityconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetEqualityConstraintsRoutine() has not been called");
424:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeEqualityConstraints");
425:   PetscLogEventBegin(Tao_ConstraintsEval,tao,X,CE,NULL);
426:   PetscStackPush("Tao equality constraints evaluation routine");
427:   (*tao->ops->computeequalityconstraints)(tao,X,CE,tao->user_con_equalityP);
428:   PetscStackPop;
429:   PetscLogEventEnd(Tao_ConstraintsEval,tao,X,CE,NULL);
430:   tao->nconstraints++;
431:   return(0);
432: }


435: /*@C
436:    TaoComputeInequalityConstraints - Compute the variable bounds using the
437:    routine set by TaoSetInequalityConstraintsRoutine().

439:    Collective on Tao

441:    Input Parameters:
442: .  tao - the Tao context

444:    Level: developer

446: .seealso: TaoSetInequalityConstraintsRoutine(), TaoComputeJacobianInequality()
447: @*/

449: PetscErrorCode TaoComputeInequalityConstraints(Tao tao, Vec X, Vec CI)
450: {


460:   if (!tao->ops->computeinequalityconstraints) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInequalityConstraintsRoutine() has not been called");
461:   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetInitialVector must be called before TaoComputeInequalityConstraints");
462:   PetscLogEventBegin(Tao_ConstraintsEval,tao,X,CI,NULL);
463:   PetscStackPush("Tao inequality constraints evaluation routine");
464:   (*tao->ops->computeinequalityconstraints)(tao,X,CI,tao->user_con_inequalityP);
465:   PetscStackPop;
466:   PetscLogEventEnd(Tao_ConstraintsEval,tao,X,CI,NULL);
467:   tao->nconstraints++;
468:   return(0);
469: }