Actual source code: almm.c

  1: #include <../src/tao/constrained/impls/almm/almm.h>
  2: #include <petsctao.h>
  3: #include <petsc/private/petscimpl.h>
  4: #include <petsc/private/vecimpl.h>

  6: static PetscErrorCode TaoALMMCombinePrimal_Private(Tao, Vec, Vec, Vec);
  7: static PetscErrorCode TaoALMMCombineDual_Private(Tao, Vec, Vec, Vec);
  8: static PetscErrorCode TaoALMMSplitPrimal_Private(Tao, Vec, Vec, Vec);
  9: static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao);
 10: static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao);
 11: static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao);

 13: static PetscErrorCode TaoSolve_ALMM(Tao tao)
 14: {
 15:   TAO_ALMM          *auglag = (TAO_ALMM *)tao->data;
 16:   TaoConvergedReason reason;
 17:   PetscReal          updated;

 19:   PetscFunctionBegin;
 20:   /* reset initial multiplier/slack guess */
 21:   if (!tao->recycle) {
 22:     if (tao->ineq_constrained) {
 23:       PetscCall(VecZeroEntries(auglag->Ps));
 24:       PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P));
 25:       PetscCall(VecSet(auglag->Yi, 1.0));
 26:     }
 27:     if (tao->eq_constrained) PetscCall(VecSet(auglag->Ye, 1.0));
 28:   }

 30:   /* compute initial nonlinear Lagrangian and its derivatives */
 31:   PetscCall((*auglag->sub_obj)(tao));
 32:   PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
 33:   /* print initial step and check convergence */
 34:   PetscCall(PetscInfo(tao, "Solving with %s formulation\n", TaoALMMTypes[auglag->type]));
 35:   PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
 36:   PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0));
 37:   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
 38:   /* set initial penalty factor and inner solver tolerance */
 39:   switch (auglag->type) {
 40:   case TAO_ALMM_CLASSIC:
 41:     auglag->mu = auglag->mu0;
 42:     break;
 43:   case TAO_ALMM_PHR:
 44:     auglag->cenorm = 0.0;
 45:     if (tao->eq_constrained) PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
 46:     auglag->cinorm = 0.0;
 47:     if (tao->ineq_constrained) {
 48:       PetscCall(VecCopy(auglag->Ci, auglag->Ciwork));
 49:       PetscCall(VecScale(auglag->Ciwork, -1.0));
 50:       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
 51:       PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm));
 52:     }
 53:     /* determine initial penalty factor based on the balance of constraint violation and objective function value */
 54:     auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0 * PetscAbsReal(auglag->fval) / (auglag->cenorm + auglag->cinorm)));
 55:     break;
 56:   default:
 57:     break;
 58:   }
 59:   auglag->gtol = auglag->gtol0;
 60:   PetscCall(PetscInfo(tao, "Initial penalty: %.2g\n", (double)auglag->mu));

 62:   /* start aug-lag outer loop */
 63:   while (tao->reason == TAO_CONTINUE_ITERATING) {
 64:     ++tao->niter;
 65:     /* update subsolver tolerance */
 66:     PetscCall(PetscInfo(tao, "Subsolver tolerance: ||G|| <= %e\n", (double)auglag->gtol));
 67:     PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
 68:     /* solve the bound-constrained or unconstrained subproblem */
 69:     PetscCall(VecCopy(auglag->P, auglag->subsolver->solution));
 70:     PetscCall(TaoSolve(auglag->subsolver));
 71:     PetscCall(VecCopy(auglag->subsolver->solution, auglag->P));
 72:     PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason));
 73:     tao->ksp_its += auglag->subsolver->ksp_its;
 74:     if (reason != TAO_CONVERGED_GATOL) PetscCall(PetscInfo(tao, "Subsolver failed to converge, reason: %s\n", TaoConvergedReasons[reason]));
 75:     /* evaluate solution and test convergence */
 76:     PetscCall((*auglag->sub_obj)(tao));
 77:     PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
 78:     /* decide whether to update multipliers or not */
 79:     updated = 0.0;
 80:     if (auglag->cnorm <= auglag->ytol) {
 81:       PetscCall(PetscInfo(tao, "Multipliers updated: ||C|| <= %e\n", (double)auglag->ytol));
 82:       /* constraints are good, update multipliers and convergence tolerances */
 83:       if (tao->eq_constrained) {
 84:         PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
 85:         PetscCall(VecSet(auglag->Cework, auglag->ye_max));
 86:         PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
 87:         PetscCall(VecSet(auglag->Cework, auglag->ye_min));
 88:         PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
 89:       }
 90:       if (tao->ineq_constrained) {
 91:         PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
 92:         PetscCall(VecSet(auglag->Ciwork, auglag->yi_max));
 93:         PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
 94:         PetscCall(VecSet(auglag->Ciwork, auglag->yi_min));
 95:         PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi));
 96:       }
 97:       /* tolerances are updated only for non-PHR methods */
 98:       if (auglag->type != TAO_ALMM_PHR) {
 99:         auglag->ytol = PetscMax(tao->catol, auglag->ytol / PetscPowReal(auglag->mu, auglag->mu_pow_good));
100:         auglag->gtol = PetscMax(tao->gatol, auglag->gtol / auglag->mu);
101:       }
102:       updated = 1.0;
103:     } else {
104:       /* constraints are bad, update penalty factor */
105:       auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac * auglag->mu);
106:       /* tolerances are reset only for non-PHR methods */
107:       if (auglag->type != TAO_ALMM_PHR) {
108:         auglag->ytol = PetscMax(tao->catol, 0.1 / PetscPowReal(auglag->mu, auglag->mu_pow_bad));
109:         auglag->gtol = PetscMax(tao->gatol, 1.0 / auglag->mu);
110:       }
111:       PetscCall(PetscInfo(tao, "Penalty increased: mu = %.2g\n", (double)auglag->mu));
112:     }
113:     PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
114:     PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
115:     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
116:   }

118:   PetscFunctionReturn(PETSC_SUCCESS);
119: }

121: static PetscErrorCode TaoView_ALMM(Tao tao, PetscViewer viewer)
122: {
123:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
124:   PetscBool isascii;

126:   PetscFunctionBegin;
127:   PetscCall(PetscViewerASCIIPushTab(viewer));
128:   PetscCall(TaoView(auglag->subsolver, viewer));
129:   PetscCall(PetscViewerASCIIPopTab(viewer));
130:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
131:   if (isascii) {
132:     PetscCall(PetscViewerASCIIPushTab(viewer));
133:     PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
134:     PetscCall(PetscViewerASCIIPopTab(viewer));
135:   }
136:   PetscFunctionReturn(PETSC_SUCCESS);
137: }

139: static PetscErrorCode TaoSetUp_ALMM(Tao tao)
140: {
141:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
142:   VecType   vec_type;
143:   Vec       SL, SU;
144:   PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE;

146:   PetscFunctionBegin;
147:   PetscCheck(!tao->ineq_doublesided, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TAOALMM does not support double-sided inequality constraint definition. Please restructure your inequality constraint to fit the form c(x) >= 0.");
148:   PetscCheck(tao->eq_constrained || tao->ineq_constrained, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "Equality and/or inequality constraints must be defined before solver setup.");
149:   PetscCall(TaoComputeVariableBounds(tao));
150:   /* alias base vectors and create extras */
151:   PetscCall(VecGetType(tao->solution, &vec_type));
152:   auglag->Px = tao->solution;
153:   if (!tao->gradient) { /* base gradient */
154:     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
155:   }
156:   auglag->LgradX = tao->gradient;
157:   if (!auglag->Xwork) { /* opt var work vector */
158:     PetscCall(VecDuplicate(tao->solution, &auglag->Xwork));
159:   }
160:   if (tao->eq_constrained) {
161:     auglag->Ce = tao->constraints_equality;
162:     auglag->Ae = tao->jacobian_equality;
163:     if (!auglag->Ye) { /* equality multipliers */
164:       PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye));
165:     }
166:     if (!auglag->Cework) PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework));
167:   }
168:   if (tao->ineq_constrained) {
169:     auglag->Ci = tao->constraints_inequality;
170:     auglag->Ai = tao->jacobian_inequality;
171:     if (!auglag->Yi) { /* inequality multipliers */
172:       PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi));
173:     }
174:     if (!auglag->Ciwork) PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork));
175:     if (!auglag->Cizero) {
176:       PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero));
177:       PetscCall(VecZeroEntries(auglag->Cizero));
178:     }
179:     if (!auglag->Ps) { /* slack vars */
180:       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps));
181:     }
182:     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
183:       PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS));
184:     }
185:     /* create vector for combined primal space and the associated communication objects */
186:     if (!auglag->P) {
187:       PetscCall(PetscMalloc1(2, &auglag->Parr));
188:       auglag->Parr[0] = auglag->Px;
189:       auglag->Parr[1] = auglag->Ps;
190:       PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
191:       PetscCall(PetscMalloc1(2, &auglag->Pscatter));
192:       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
193:       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]));
194:     }
195:     if (tao->eq_constrained) {
196:       /* create vector for combined dual space and the associated communication objects */
197:       if (!auglag->Y) {
198:         PetscCall(PetscMalloc1(2, &auglag->Yarr));
199:         auglag->Yarr[0] = auglag->Ye;
200:         auglag->Yarr[1] = auglag->Yi;
201:         PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
202:         PetscCall(PetscMalloc1(2, &auglag->Yscatter));
203:         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
204:         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
205:       }
206:       if (!auglag->C) PetscCall(VecDuplicate(auglag->Y, &auglag->C));
207:     } else {
208:       if (!auglag->C) auglag->C = auglag->Ci;
209:       if (!auglag->Y) auglag->Y = auglag->Yi;
210:     }
211:   } else {
212:     if (!auglag->P) auglag->P = auglag->Px;
213:     if (!auglag->G) auglag->G = auglag->LgradX;
214:     if (!auglag->C) auglag->C = auglag->Ce;
215:     if (!auglag->Y) auglag->Y = auglag->Ye;
216:   }
217:   /* create tao primal solution and gradient to interface with subsolver */
218:   PetscCall(VecDuplicate(auglag->P, &auglag->Psub));
219:   PetscCall(VecDuplicate(auglag->P, &auglag->G));
220:   /* initialize parameters */
221:   if (auglag->type == TAO_ALMM_PHR) {
222:     auglag->mu_fac = 10.0;
223:     auglag->yi_min = 0.0;
224:     auglag->ytol0  = 0.5;
225:     auglag->gtol0  = tao->gatol;
226:     if (tao->gatol_changed && tao->catol_changed) {
227:       PetscCall(PetscInfo(tao, "TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n"));
228:       tao->catol = tao->gatol;
229:     }
230:   }
231:   /* set the Lagrangian formulation type for the subsolver */
232:   switch (auglag->type) {
233:   case TAO_ALMM_CLASSIC:
234:     auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
235:     break;
236:   case TAO_ALMM_PHR:
237:     auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private;
238:     break;
239:   default:
240:     break;
241:   }
242:   /* set up the subsolver */
243:   PetscCall(TaoSetSolution(auglag->subsolver, auglag->Psub));
244:   PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag));
245:   PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag));
246:   if (tao->bounded) {
247:     /* make sure that the subsolver is a bound-constrained method */
248:     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
249:     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
250:     if (is_cg) {
251:       PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
252:       PetscCall(PetscInfo(tao, "TAOCG detected for bound-constrained problem, switching to TAOBNCG.\n"));
253:     }
254:     if (is_lmvm) {
255:       PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
256:       PetscCall(PetscInfo(tao, "TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS.\n"));
257:     }
258:     /* create lower and upper bound clone vectors for subsolver */
259:     if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL));
260:     if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU));
261:     if (tao->ineq_constrained) {
262:       /* create lower and upper bounds for slack, set lower to 0 */
263:       PetscCall(VecDuplicate(auglag->Ci, &SL));
264:       PetscCall(VecSet(SL, 0.0));
265:       PetscCall(VecDuplicate(auglag->Ci, &SU));
266:       PetscCall(VecSet(SU, PETSC_INFINITY));
267:       /* combine opt var bounds with slack bounds */
268:       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
269:       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
270:       /* destroy work vectors */
271:       PetscCall(VecDestroy(&SL));
272:       PetscCall(VecDestroy(&SU));
273:     } else {
274:       /* no inequality constraints, just copy bounds into the subsolver */
275:       PetscCall(VecCopy(tao->XL, auglag->PL));
276:       PetscCall(VecCopy(tao->XU, auglag->PU));
277:     }
278:     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
279:   }
280:   PetscCall(TaoSetUp(auglag->subsolver));

282:   PetscFunctionReturn(PETSC_SUCCESS);
283: }

285: static PetscErrorCode TaoDestroy_ALMM(Tao tao)
286: {
287:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

289:   PetscFunctionBegin;
290:   PetscCall(TaoDestroy(&auglag->subsolver));
291:   PetscCall(VecDestroy(&auglag->Psub));
292:   PetscCall(VecDestroy(&auglag->G));
293:   if (tao->setupcalled) {
294:     PetscCall(VecDestroy(&auglag->Xwork));
295:     if (tao->eq_constrained) {
296:       PetscCall(VecDestroy(&auglag->Ye));     /* equality multipliers */
297:       PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */
298:     }
299:     if (tao->ineq_constrained) {
300:       PetscCall(VecDestroy(&auglag->Ps));     /* slack vars */
301:       PetscCall(PetscFree(auglag->Parr));     /* array of primal vectors */
302:       PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */
303:       PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */
304:       PetscCall(VecDestroy(&auglag->Yi));     /* inequality multipliers */
305:       PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */
306:       PetscCall(VecDestroy(&auglag->P));      /* combo primal solution */
307:       PetscCall(ISDestroy(&auglag->Pis[0]));  /* index set for X inside P */
308:       PetscCall(ISDestroy(&auglag->Pis[1]));  /* index set for S inside P */
309:       PetscCall(PetscFree(auglag->Pis));      /* array of P index sets */
310:       PetscCall(VecScatterDestroy(&auglag->Pscatter[0]));
311:       PetscCall(VecScatterDestroy(&auglag->Pscatter[1]));
312:       PetscCall(PetscFree(auglag->Pscatter));
313:       if (tao->eq_constrained) {
314:         PetscCall(VecDestroy(&auglag->Y));     /* combo multipliers */
315:         PetscCall(PetscFree(auglag->Yarr));    /* array of dual vectors */
316:         PetscCall(VecDestroy(&auglag->C));     /* combo constraints */
317:         PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */
318:         PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */
319:         PetscCall(PetscFree(auglag->Yis));
320:         PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
321:         PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
322:         PetscCall(PetscFree(auglag->Yscatter));
323:       }
324:     }
325:     if (tao->bounded) {
326:       PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */
327:       PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */
328:     }
329:   }
330:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", NULL));
331:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", NULL));
332:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", NULL));
333:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", NULL));
334:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", NULL));
335:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", NULL));
336:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", NULL));
337:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", NULL));
338:   PetscCall(PetscFree(tao->data));
339:   PetscFunctionReturn(PETSC_SUCCESS);
340: }

342: static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao, PetscOptionItems *PetscOptionsObject)
343: {
344:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
345:   PetscInt  i;

347:   PetscFunctionBegin;
348:   PetscOptionsHeadBegin(PetscOptionsObject, "Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems.");
349:   PetscCall(PetscOptionsReal("-tao_almm_mu_init", "initial penalty parameter", "", auglag->mu0, &auglag->mu0, NULL));
350:   PetscCall(PetscOptionsReal("-tao_almm_mu_factor", "increase factor for the penalty parameter", "", auglag->mu_fac, &auglag->mu_fac, NULL));
351:   PetscCall(PetscOptionsReal("-tao_almm_mu_power_good", "exponential for penalty parameter when multiplier update is accepted", "", auglag->mu_pow_good, &auglag->mu_pow_good, NULL));
352:   PetscCall(PetscOptionsReal("-tao_almm_mu_power_bad", "exponential for penalty parameter when multiplier update is rejected", "", auglag->mu_pow_bad, &auglag->mu_pow_bad, NULL));
353:   PetscCall(PetscOptionsReal("-tao_almm_mu_max", "maximum safeguard for penalty parameter updates", "", auglag->mu_max, &auglag->mu_max, NULL));
354:   PetscCall(PetscOptionsReal("-tao_almm_ye_min", "minimum safeguard for equality multiplier updates", "", auglag->ye_min, &auglag->ye_min, NULL));
355:   PetscCall(PetscOptionsReal("-tao_almm_ye_max", "maximum safeguard for equality multipliers updates", "", auglag->ye_max, &auglag->ye_max, NULL));
356:   PetscCall(PetscOptionsReal("-tao_almm_yi_min", "minimum safeguard for inequality multipliers updates", "", auglag->yi_min, &auglag->yi_min, NULL));
357:   PetscCall(PetscOptionsReal("-tao_almm_yi_max", "maximum safeguard for inequality multipliers updates", "", auglag->yi_max, &auglag->yi_max, NULL));
358:   PetscCall(PetscOptionsEnum("-tao_almm_type", "augmented Lagrangian formulation type for the subproblem", "TaoALMMType", TaoALMMTypes, (PetscEnum)auglag->type, (PetscEnum *)&auglag->type, NULL));
359:   PetscOptionsHeadEnd();
360:   PetscCall(TaoSetOptionsPrefix(auglag->subsolver, ((PetscObject)tao)->prefix));
361:   PetscCall(TaoAppendOptionsPrefix(auglag->subsolver, "tao_almm_subsolver_"));
362:   PetscCall(TaoSetFromOptions(auglag->subsolver));
363:   for (i = 0; i < tao->numbermonitors; i++) {
364:     PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
365:     PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
366:     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) auglag->info = PETSC_TRUE;
367:   }
368:   PetscFunctionReturn(PETSC_SUCCESS);
369: }

371: /* -------------------------------------------------------- */

373: /*MC
374:   TAOALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.

376:   Options Database Keys:
377: + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
378: . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
379: . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
380: . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
381: . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
382: . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
383: . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
384: . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
385: . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
386: - -tao_almm_type <phr,classic>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: phr)

388:   Level: beginner

390:   Notes:
391:   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
392:   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.

394:   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
395:   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
396:   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation may converge
397:   faster for smaller problems but is highly susceptible to poor step lengths in the subproblem due to the positivity
398:   constraint on slack variables. PHR avoids this issue by eliminating the slack variables entirely, and is highly
399:   desirable for problems with a large number of inequality constraints.

401:   The subproblem is solved using a nested first-order TAO solver (default: `TAOBQNLS`). The user can retrieve a
402:   pointer to the subsolver via `TaoALMMGetSubsolver()` or pass command line arguments to it using the
403:   "-tao_almm_subsolver_" prefix. Currently, `TAOALMM` does not support second-order methods for the subproblem.

405: .vb
406:   while unconverged
407:     solve argmin_x L(x) s.t. l <= x <= u
408:     if ||c|| <= y_tol
409:       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
410:         problem converged, return solution
411:       else
412:         constraints sufficiently improved
413:         update multipliers and tighten tolerances
414:       endif
415:     else
416:       constraints did not improve
417:       update penalty and loosen tolerances
418:     endif
419:   endwhile
420: .ve

422: .seealso: `TAOALMM`, `Tao`, `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`,
423:           `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()`
424: M*/
425: PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao)
426: {
427:   TAO_ALMM *auglag;

429:   PetscFunctionBegin;
430:   PetscCall(PetscNew(&auglag));

432:   tao->ops->destroy        = TaoDestroy_ALMM;
433:   tao->ops->setup          = TaoSetUp_ALMM;
434:   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
435:   tao->ops->view           = TaoView_ALMM;
436:   tao->ops->solve          = TaoSolve_ALMM;

438:   tao->gatol = 1.e-5;
439:   tao->grtol = 0.0;
440:   tao->gttol = 0.0;
441:   tao->catol = 1.e-5;
442:   tao->crtol = 0.0;

444:   tao->data           = (void *)auglag;
445:   auglag->parent      = tao;
446:   auglag->mu0         = 10.0;
447:   auglag->mu          = auglag->mu0;
448:   auglag->mu_fac      = 10.0;
449:   auglag->mu_max      = PETSC_INFINITY;
450:   auglag->mu_pow_good = 0.9;
451:   auglag->mu_pow_bad  = 0.1;
452:   auglag->ye_min      = PETSC_NINFINITY;
453:   auglag->ye_max      = PETSC_INFINITY;
454:   auglag->yi_min      = PETSC_NINFINITY;
455:   auglag->yi_max      = PETSC_INFINITY;
456:   auglag->ytol0       = 0.1 / PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
457:   auglag->ytol        = auglag->ytol0;
458:   auglag->gtol0       = 1.0 / auglag->mu0;
459:   auglag->gtol        = auglag->gtol0;

461:   auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
462:   auglag->type    = TAO_ALMM_PHR;
463:   auglag->info    = PETSC_FALSE;

465:   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &auglag->subsolver));
466:   PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
467:   PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
468:   PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000));
469:   PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
470:   PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
471:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver, (PetscObject)tao, 1));

473:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
474:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
475:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
476:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
477:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
478:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
479:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
480:   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
481:   PetscFunctionReturn(PETSC_SUCCESS);
482: }

484: static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P)
485: {
486:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

488:   PetscFunctionBegin;
489:   if (tao->ineq_constrained) {
490:     PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
491:     PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
492:     PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
493:     PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
494:   } else {
495:     PetscCall(VecCopy(X, P));
496:   }
497:   PetscFunctionReturn(PETSC_SUCCESS);
498: }

500: static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y)
501: {
502:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

504:   PetscFunctionBegin;
505:   if (tao->eq_constrained) {
506:     if (tao->ineq_constrained) {
507:       PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
508:       PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
509:       PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
510:       PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
511:     } else {
512:       PetscCall(VecCopy(EQ, Y));
513:     }
514:   } else {
515:     PetscCall(VecCopy(IN, Y));
516:   }
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S)
521: {
522:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

524:   PetscFunctionBegin;
525:   if (tao->ineq_constrained) {
526:     PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
527:     PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
528:     PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
529:     PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
530:   } else {
531:     PetscCall(VecCopy(P, X));
532:   }
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }

536: /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
537: static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao)
538: {
539:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

541:   PetscFunctionBegin;
542:   /* if bounded, project the gradient */
543:   if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
544:   if (auglag->type == TAO_ALMM_PHR) {
545:     PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
546:     auglag->cenorm = 0.0;
547:     if (tao->eq_constrained) PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
548:     auglag->cinorm = 0.0;
549:     if (tao->ineq_constrained) {
550:       PetscCall(VecCopy(auglag->Yi, auglag->Ciwork));
551:       PetscCall(VecScale(auglag->Ciwork, -1.0 / auglag->mu));
552:       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
553:       PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
554:     }
555:     auglag->cnorm_old = auglag->cnorm;
556:     auglag->cnorm     = PetscMax(auglag->cenorm, auglag->cinorm);
557:     auglag->ytol      = auglag->ytol0 * auglag->cnorm_old;
558:   } else {
559:     PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
560:     PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
561:   }
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao)
566: {
567:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;

569:   PetscFunctionBegin;
570:   /* split solution into primal and slack components */
571:   PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));

573:   /* compute f, df/dx and the constraints */
574:   PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
575:   if (tao->eq_constrained) {
576:     PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
577:     PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
578:   }
579:   if (tao->ineq_constrained) {
580:     PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
581:     PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
582:     switch (auglag->type) {
583:     case TAO_ALMM_CLASSIC:
584:       /* classic formulation converts inequality to equality constraints via slack variables */
585:       PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
586:       break;
587:     case TAO_ALMM_PHR:
588:       /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
589:       PetscCall(VecScale(auglag->Ci, -1.0));
590:       PetscCall(MatScale(auglag->Ai, -1.0));
591:       break;
592:     default:
593:       break;
594:     }
595:   }
596:   /* combine constraints into one vector */
597:   PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
598:   PetscFunctionReturn(PETSC_SUCCESS);
599: }

601: /*
602: Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)]

604: dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai]

606: dLphr/dS = 0
607: */
608: static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao)
609: {
610:   TAO_ALMM *auglag  = (TAO_ALMM *)tao->data;
611:   PetscReal eq_norm = 0.0, ineq_norm = 0.0;

613:   PetscFunctionBegin;
614:   PetscCall(TaoALMMEvaluateIterate_Private(tao));
615:   if (tao->eq_constrained) {
616:     /* Ce_work = mu*(Ce + Ye/mu) */
617:     PetscCall(VecWAXPY(auglag->Cework, 1.0 / auglag->mu, auglag->Ye, auglag->Ce));
618:     PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
619:     PetscCall(VecScale(auglag->Cework, auglag->mu));
620:     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
621:     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
622:   }
623:   if (tao->ineq_constrained) {
624:     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
625:     PetscCall(VecWAXPY(auglag->Ciwork, 1.0 / auglag->mu, auglag->Yi, auglag->Ci));
626:     PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
627:     PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
628:     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
629:     PetscCall(VecScale(auglag->Ciwork, auglag->mu));
630:     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
631:     /* dL/dS = 0 because there are no slacks in PHR */
632:     PetscCall(VecZeroEntries(auglag->LgradS));
633:   }
634:   /* combine gradient together */
635:   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
636:   /* compute L = f + 0.5 * mu * [(Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)] */
637:   auglag->Lval = auglag->fval + 0.5 * auglag->mu * (eq_norm + ineq_norm);
638:   PetscFunctionReturn(PETSC_SUCCESS);
639: }

641: /*
642: Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]

644: dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi]

646: dLc/dS = -[Yi + mu*(Ci - S)]
647: */
648: static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao)
649: {
650:   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
651:   PetscReal yeTce = 0.0, yiTcims = 0.0, ceTce = 0.0, cimsTcims = 0.0;

653:   PetscFunctionBegin;
654:   PetscCall(TaoALMMEvaluateIterate_Private(tao));
655:   if (tao->eq_constrained) {
656:     /* compute scalar contributions */
657:     PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce));
658:     PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce));
659:     /* dL/dX += ye^T Ae */
660:     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
661:     /* dL/dX += 0.5 * mu * ce^T Ae */
662:     PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
663:     PetscCall(VecAXPY(auglag->LgradX, 0.5 * auglag->mu, auglag->Xwork));
664:   }
665:   if (tao->ineq_constrained) {
666:     /* compute scalar contributions */
667:     PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
668:     PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
669:     /* dL/dX += yi^T Ai */
670:     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
671:     /* dL/dX += 0.5 * mu * (ci - s)^T Ai */
672:     PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
673:     PetscCall(VecAXPY(auglag->LgradX, 0.5 * auglag->mu, auglag->Xwork));
674:     /* dL/dS = -[yi + mu*(ci - s)] */
675:     PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
676:     PetscCall(VecScale(auglag->LgradS, -1.0));
677:   }
678:   /* combine gradient together */
679:   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
680:   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
681:   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5 * auglag->mu * (ceTce + cimsTcims);
682:   PetscFunctionReturn(PETSC_SUCCESS);
683: }

685: PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx)
686: {
687:   TAO_ALMM *auglag = (TAO_ALMM *)ctx;

689:   PetscFunctionBegin;
690:   PetscCall(VecCopy(P, auglag->P));
691:   PetscCall((*auglag->sub_obj)(auglag->parent));
692:   *Lval = auglag->Lval;
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: }

696: PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx)
697: {
698:   TAO_ALMM *auglag = (TAO_ALMM *)ctx;

700:   PetscFunctionBegin;
701:   PetscCall(VecCopy(P, auglag->P));
702:   PetscCall((*auglag->sub_obj)(auglag->parent));
703:   PetscCall(VecCopy(auglag->G, G));
704:   *Lval = auglag->Lval;
705:   PetscFunctionReturn(PETSC_SUCCESS);
706: }